Merge branch 'origin/release-2020' into merge-release-2020-into-master
[alexxy/gromacs.git] / src / gromacs / mdlib / sim_util.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013-2019,2020, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #include <cmath>
42 #include <cstdint>
43 #include <cstdio>
44 #include <cstring>
45
46 #include <array>
47 #include <optional>
48
49 #include "gromacs/applied_forces/awh/awh.h"
50 #include "gromacs/domdec/dlbtiming.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/domdec/gpuhaloexchange.h"
54 #include "gromacs/domdec/partition.h"
55 #include "gromacs/essentialdynamics/edsam.h"
56 #include "gromacs/ewald/pme.h"
57 #include "gromacs/ewald/pme_pp.h"
58 #include "gromacs/ewald/pme_pp_comm_gpu.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
61 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
62 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
63 #include "gromacs/gpu_utils/gpu_utils.h"
64 #include "gromacs/imd/imd.h"
65 #include "gromacs/listed_forces/disre.h"
66 #include "gromacs/listed_forces/gpubonded.h"
67 #include "gromacs/listed_forces/listed_forces.h"
68 #include "gromacs/listed_forces/orires.h"
69 #include "gromacs/math/arrayrefwithpadding.h"
70 #include "gromacs/math/functions.h"
71 #include "gromacs/math/units.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/math/vecdump.h"
74 #include "gromacs/mdlib/calcmu.h"
75 #include "gromacs/mdlib/calcvir.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/dispersioncorrection.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/force.h"
80 #include "gromacs/mdlib/force_flags.h"
81 #include "gromacs/mdlib/forcerec.h"
82 #include "gromacs/mdlib/gmx_omp_nthreads.h"
83 #include "gromacs/mdlib/update.h"
84 #include "gromacs/mdlib/vsite.h"
85 #include "gromacs/mdlib/wall.h"
86 #include "gromacs/mdlib/wholemoleculetransform.h"
87 #include "gromacs/mdtypes/commrec.h"
88 #include "gromacs/mdtypes/enerdata.h"
89 #include "gromacs/mdtypes/forcebuffers.h"
90 #include "gromacs/mdtypes/forceoutput.h"
91 #include "gromacs/mdtypes/forcerec.h"
92 #include "gromacs/mdtypes/iforceprovider.h"
93 #include "gromacs/mdtypes/inputrec.h"
94 #include "gromacs/mdtypes/md_enums.h"
95 #include "gromacs/mdtypes/mdatom.h"
96 #include "gromacs/mdtypes/multipletimestepping.h"
97 #include "gromacs/mdtypes/simulation_workload.h"
98 #include "gromacs/mdtypes/state.h"
99 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
100 #include "gromacs/nbnxm/gpu_data_mgmt.h"
101 #include "gromacs/nbnxm/nbnxm.h"
102 #include "gromacs/nbnxm/nbnxm_gpu.h"
103 #include "gromacs/pbcutil/ishift.h"
104 #include "gromacs/pbcutil/pbc.h"
105 #include "gromacs/pulling/pull.h"
106 #include "gromacs/pulling/pull_rotation.h"
107 #include "gromacs/timing/cyclecounter.h"
108 #include "gromacs/timing/gpu_timing.h"
109 #include "gromacs/timing/wallcycle.h"
110 #include "gromacs/timing/wallcyclereporting.h"
111 #include "gromacs/timing/walltime_accounting.h"
112 #include "gromacs/topology/topology.h"
113 #include "gromacs/utility/arrayref.h"
114 #include "gromacs/utility/basedefinitions.h"
115 #include "gromacs/utility/cstringutil.h"
116 #include "gromacs/utility/exceptions.h"
117 #include "gromacs/utility/fatalerror.h"
118 #include "gromacs/utility/fixedcapacityvector.h"
119 #include "gromacs/utility/gmxassert.h"
120 #include "gromacs/utility/gmxmpi.h"
121 #include "gromacs/utility/logger.h"
122 #include "gromacs/utility/smalloc.h"
123 #include "gromacs/utility/strconvert.h"
124 #include "gromacs/utility/sysinfo.h"
125
126 #include "gpuforcereduction.h"
127
128 using gmx::ArrayRef;
129 using gmx::AtomLocality;
130 using gmx::DomainLifetimeWorkload;
131 using gmx::ForceOutputs;
132 using gmx::ForceWithShiftForces;
133 using gmx::InteractionLocality;
134 using gmx::RVec;
135 using gmx::SimulationWorkload;
136 using gmx::StepWorkload;
137
138 // TODO: this environment variable allows us to verify before release
139 // that on less common architectures the total cost of polling is not larger than
140 // a blocking wait (so polling does not introduce overhead when the static
141 // PME-first ordering would suffice).
142 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
143
144 static void sum_forces(ArrayRef<RVec> f, ArrayRef<const RVec> forceToAdd)
145 {
146     GMX_ASSERT(f.size() >= forceToAdd.size(), "Accumulation buffer should be sufficiently large");
147     const int end = forceToAdd.size();
148
149     int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
150 #pragma omp parallel for num_threads(nt) schedule(static)
151     for (int i = 0; i < end; i++)
152     {
153         rvec_inc(f[i], forceToAdd[i]);
154     }
155 }
156
157 static void calc_virial(int                              start,
158                         int                              homenr,
159                         const rvec                       x[],
160                         const gmx::ForceWithShiftForces& forceWithShiftForces,
161                         tensor                           vir_part,
162                         const matrix                     box,
163                         t_nrnb*                          nrnb,
164                         const t_forcerec*                fr,
165                         PbcType                          pbcType)
166 {
167     /* The short-range virial from surrounding boxes */
168     const rvec* fshift = as_rvec_array(forceWithShiftForces.shiftForces().data());
169     calc_vir(SHIFTS, fr->shift_vec, fshift, vir_part, pbcType == PbcType::Screw, box);
170     inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
171
172     /* Calculate partial virial, for local atoms only, based on short range.
173      * Total virial is computed in global_stat, called from do_md
174      */
175     const rvec* f = as_rvec_array(forceWithShiftForces.force().data());
176     f_calc_vir(start, start + homenr, x, f, vir_part, box);
177     inc_nrnb(nrnb, eNR_VIRIAL, homenr);
178
179     if (debug)
180     {
181         pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
182     }
183 }
184
185 static void pull_potential_wrapper(const t_commrec*               cr,
186                                    const t_inputrec*              ir,
187                                    const matrix                   box,
188                                    gmx::ArrayRef<const gmx::RVec> x,
189                                    gmx::ForceWithVirial*          force,
190                                    const t_mdatoms*               mdatoms,
191                                    gmx_enerdata_t*                enerd,
192                                    pull_t*                        pull_work,
193                                    const real*                    lambda,
194                                    double                         t,
195                                    gmx_wallcycle_t                wcycle)
196 {
197     t_pbc pbc;
198     real  dvdl;
199
200     /* Calculate the center of mass forces, this requires communication,
201      * which is why pull_potential is called close to other communication.
202      */
203     wallcycle_start(wcycle, ewcPULLPOT);
204     set_pbc(&pbc, ir->pbcType, box);
205     dvdl = 0;
206     enerd->term[F_COM_PULL] +=
207             pull_potential(pull_work, mdatoms->massT, &pbc, cr, t, lambda[efptRESTRAINT],
208                            as_rvec_array(x.data()), force, &dvdl);
209     enerd->dvdl_lin[efptRESTRAINT] += dvdl;
210     wallcycle_stop(wcycle, ewcPULLPOT);
211 }
212
213 static void pme_receive_force_ener(t_forcerec*           fr,
214                                    const t_commrec*      cr,
215                                    gmx::ForceWithVirial* forceWithVirial,
216                                    gmx_enerdata_t*       enerd,
217                                    bool                  useGpuPmePpComms,
218                                    bool                  receivePmeForceToGpu,
219                                    gmx_wallcycle_t       wcycle)
220 {
221     real  e_q, e_lj, dvdl_q, dvdl_lj;
222     float cycles_ppdpme, cycles_seppme;
223
224     cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
225     dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
226
227     /* In case of node-splitting, the PP nodes receive the long-range
228      * forces, virial and energy from the PME nodes here.
229      */
230     wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
231     dvdl_q  = 0;
232     dvdl_lj = 0;
233     gmx_pme_receive_f(fr->pmePpCommGpu.get(), cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
234                       useGpuPmePpComms, receivePmeForceToGpu, &cycles_seppme);
235     enerd->term[F_COUL_RECIP] += e_q;
236     enerd->term[F_LJ_RECIP] += e_lj;
237     enerd->dvdl_lin[efptCOUL] += dvdl_q;
238     enerd->dvdl_lin[efptVDW] += dvdl_lj;
239
240     if (wcycle)
241     {
242         dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
243     }
244     wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
245 }
246
247 static void print_large_forces(FILE*                fp,
248                                const t_mdatoms*     md,
249                                const t_commrec*     cr,
250                                int64_t              step,
251                                real                 forceTolerance,
252                                ArrayRef<const RVec> x,
253                                ArrayRef<const RVec> f)
254 {
255     real       force2Tolerance = gmx::square(forceTolerance);
256     gmx::index numNonFinite    = 0;
257     for (int i = 0; i < md->homenr; i++)
258     {
259         real force2    = norm2(f[i]);
260         bool nonFinite = !std::isfinite(force2);
261         if (force2 >= force2Tolerance || nonFinite)
262         {
263             fprintf(fp, "step %" PRId64 " atom %6d  x %8.3f %8.3f %8.3f  force %12.5e\n", step,
264                     ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
265         }
266         if (nonFinite)
267         {
268             numNonFinite++;
269         }
270     }
271     if (numNonFinite > 0)
272     {
273         /* Note that with MPI this fatal call on one rank might interrupt
274          * the printing on other ranks. But we can only avoid that with
275          * an expensive MPI barrier that we would need at each step.
276          */
277         gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
278     }
279 }
280
281 //! When necessary, spreads forces on vsites and computes the virial for \p forceOutputs->forceWithShiftForces()
282 static void postProcessForceWithShiftForces(t_nrnb*                   nrnb,
283                                             gmx_wallcycle_t           wcycle,
284                                             const matrix              box,
285                                             ArrayRef<const RVec>      x,
286                                             ForceOutputs*             forceOutputs,
287                                             tensor                    vir_force,
288                                             const t_mdatoms&          mdatoms,
289                                             const t_forcerec&         fr,
290                                             gmx::VirtualSitesHandler* vsite,
291                                             const StepWorkload&       stepWork)
292 {
293     ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
294
295     /* If we have NoVirSum forces, but we do not calculate the virial,
296      * we later sum the forceWithShiftForces buffer together with
297      * the noVirSum buffer and spread the combined vsite forces at once.
298      */
299     if (vsite && (!forceOutputs->haveForceWithVirial() || stepWork.computeVirial))
300     {
301         using VirialHandling = gmx::VirtualSitesHandler::VirialHandling;
302
303         auto                 f      = forceWithShiftForces.force();
304         auto                 fshift = forceWithShiftForces.shiftForces();
305         const VirialHandling virialHandling =
306                 (stepWork.computeVirial ? VirialHandling::Pbc : VirialHandling::None);
307         vsite->spreadForces(x, f, virialHandling, fshift, nullptr, nrnb, box, wcycle);
308         forceWithShiftForces.haveSpreadVsiteForces() = true;
309     }
310
311     if (stepWork.computeVirial)
312     {
313         /* Calculation of the virial must be done after vsites! */
314         calc_virial(0, mdatoms.homenr, as_rvec_array(x.data()), forceWithShiftForces, vir_force,
315                     box, nrnb, &fr, fr.pbcType);
316     }
317 }
318
319 //! Spread, compute virial for and sum forces, when necessary
320 static void postProcessForces(const t_commrec*          cr,
321                               int64_t                   step,
322                               t_nrnb*                   nrnb,
323                               gmx_wallcycle_t           wcycle,
324                               const matrix              box,
325                               ArrayRef<const RVec>      x,
326                               ForceOutputs*             forceOutputs,
327                               tensor                    vir_force,
328                               const t_mdatoms*          mdatoms,
329                               const t_forcerec*         fr,
330                               gmx::VirtualSitesHandler* vsite,
331                               const StepWorkload&       stepWork)
332 {
333     // Extract the final output force buffer, which is also the buffer for forces with shift forces
334     ArrayRef<RVec> f = forceOutputs->forceWithShiftForces().force();
335
336     if (forceOutputs->haveForceWithVirial())
337     {
338         auto& forceWithVirial = forceOutputs->forceWithVirial();
339
340         if (vsite)
341         {
342             /* Spread the mesh force on virtual sites to the other particles...
343              * This is parallellized. MPI communication is performed
344              * if the constructing atoms aren't local.
345              */
346             GMX_ASSERT(!stepWork.computeVirial || f.data() != forceWithVirial.force_.data(),
347                        "We need separate force buffers for shift and virial forces when "
348                        "computing the virial");
349             GMX_ASSERT(!stepWork.computeVirial
350                                || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
351                        "We should spread the force with shift forces separately when computing "
352                        "the virial");
353             const gmx::VirtualSitesHandler::VirialHandling virialHandling =
354                     (stepWork.computeVirial ? gmx::VirtualSitesHandler::VirialHandling::NonLinear
355                                             : gmx::VirtualSitesHandler::VirialHandling::None);
356             matrix virial = { { 0 } };
357             vsite->spreadForces(x, forceWithVirial.force_, virialHandling, {}, virial, nrnb, box, wcycle);
358             forceWithVirial.addVirialContribution(virial);
359         }
360
361         if (stepWork.computeVirial)
362         {
363             /* Now add the forces, this is local */
364             sum_forces(f, forceWithVirial.force_);
365
366             /* Add the direct virial contributions */
367             GMX_ASSERT(
368                     forceWithVirial.computeVirial_,
369                     "forceWithVirial should request virial computation when we request the virial");
370             m_add(vir_force, forceWithVirial.getVirial(), vir_force);
371
372             if (debug)
373             {
374                 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
375             }
376         }
377     }
378     else
379     {
380         GMX_ASSERT(vsite == nullptr || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
381                    "We should have spread the vsite forces (earlier)");
382     }
383
384     if (fr->print_force >= 0)
385     {
386         print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
387     }
388 }
389
390 static void do_nb_verlet(t_forcerec*                fr,
391                          const interaction_const_t* ic,
392                          gmx_enerdata_t*            enerd,
393                          const StepWorkload&        stepWork,
394                          const InteractionLocality  ilocality,
395                          const int                  clearF,
396                          const int64_t              step,
397                          t_nrnb*                    nrnb,
398                          gmx_wallcycle_t            wcycle)
399 {
400     if (!stepWork.computeNonbondedForces)
401     {
402         /* skip non-bonded calculation */
403         return;
404     }
405
406     nonbonded_verlet_t* nbv = fr->nbv.get();
407
408     /* GPU kernel launch overhead is already timed separately */
409     if (!nbv->useGpu())
410     {
411         /* When dynamic pair-list  pruning is requested, we need to prune
412          * at nstlistPrune steps.
413          */
414         if (nbv->isDynamicPruningStepCpu(step))
415         {
416             /* Prune the pair-list beyond fr->ic->rlistPrune using
417              * the current coordinates of the atoms.
418              */
419             wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
420             nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
421             wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
422         }
423     }
424
425     nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
426 }
427
428 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
429 {
430     int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, v.ssize());
431
432     /* Note that we would like to avoid this conditional by putting it
433      * into the omp pragma instead, but then we still take the full
434      * omp parallel for overhead (at least with gcc5).
435      */
436     if (!useOpenmpThreading || nth == 1)
437     {
438         for (RVec& elem : v)
439         {
440             clear_rvec(elem);
441         }
442     }
443     else
444     {
445 #pragma omp parallel for num_threads(nth) schedule(static)
446         for (gmx::index i = 0; i < v.ssize(); i++)
447         {
448             clear_rvec(v[i]);
449         }
450     }
451 }
452
453 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
454  *
455  * \param groupOptions  Group options, containing T-coupling options
456  */
457 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
458 {
459     real nrdfCoupled   = 0;
460     real nrdfUncoupled = 0;
461     real kineticEnergy = 0;
462     for (int g = 0; g < groupOptions.ngtc; g++)
463     {
464         if (groupOptions.tau_t[g] >= 0)
465         {
466             nrdfCoupled += groupOptions.nrdf[g];
467             kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * BOLTZ;
468         }
469         else
470         {
471             nrdfUncoupled += groupOptions.nrdf[g];
472         }
473     }
474
475     /* This conditional with > also catches nrdf=0 */
476     if (nrdfCoupled > nrdfUncoupled)
477     {
478         return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
479     }
480     else
481     {
482         return 0;
483     }
484 }
485
486 /*! \brief This routine checks that the potential energy is finite.
487  *
488  * Always checks that the potential energy is finite. If step equals
489  * inputrec.init_step also checks that the magnitude of the potential energy
490  * is reasonable. Terminates with a fatal error when a check fails.
491  * Note that passing this check does not guarantee finite forces,
492  * since those use slightly different arithmetics. But in most cases
493  * there is just a narrow coordinate range where forces are not finite
494  * and energies are finite.
495  *
496  * \param[in] step      The step number, used for checking and printing
497  * \param[in] enerd     The energy data; the non-bonded group energies need to be added to
498  * enerd.term[F_EPOT] before calling this routine \param[in] inputrec  The input record
499  */
500 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
501 {
502     /* Threshold valid for comparing absolute potential energy against
503      * the kinetic energy. Normally one should not consider absolute
504      * potential energy values, but with a factor of one million
505      * we should never get false positives.
506      */
507     constexpr real c_thresholdFactor = 1e6;
508
509     bool energyIsNotFinite    = !std::isfinite(enerd.term[F_EPOT]);
510     real averageKineticEnergy = 0;
511     /* We only check for large potential energy at the initial step,
512      * because that is by far the most likely step for this too occur
513      * and because computing the average kinetic energy is not free.
514      * Note: nstcalcenergy >> 1 often does not allow to catch large energies
515      * before they become NaN.
516      */
517     if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
518     {
519         averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
520     }
521
522     if (energyIsNotFinite
523         || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
524     {
525         gmx_fatal(
526                 FARGS,
527                 "Step %" PRId64
528                 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
529                 "contributions to the energy are %g and %g, respectively. A %s potential energy "
530                 "can be caused by overlapping interactions in bonded interactions or very large%s "
531                 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
532                 "configuration, incorrect interactions or parameters in the topology.",
533                 step, enerd.term[F_EPOT], energyIsNotFinite ? "not finite" : "extremely high",
534                 enerd.term[F_LJ], enerd.term[F_COUL_SR],
535                 energyIsNotFinite ? "non-finite" : "very high", energyIsNotFinite ? " or Nan" : "");
536     }
537 }
538
539 /*! \brief Return true if there are special forces computed this step.
540  *
541  * The conditionals exactly correspond to those in computeSpecialForces().
542  */
543 static bool haveSpecialForces(const t_inputrec&          inputrec,
544                               const gmx::ForceProviders& forceProviders,
545                               const pull_t*              pull_work,
546                               const bool                 computeForces,
547                               const gmx_edsam*           ed)
548 {
549
550     return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
551             (inputrec.bPull && pull_have_potential(pull_work)) ||   // pull
552             inputrec.bRot ||                                        // enforced rotation
553             (ed != nullptr) ||                                      // flooding
554             (inputrec.bIMD && computeForces));                      // IMD
555 }
556
557 /*! \brief Compute forces and/or energies for special algorithms
558  *
559  * The intention is to collect all calls to algorithms that compute
560  * forces on local atoms only and that do not contribute to the local
561  * virial sum (but add their virial contribution separately).
562  * Eventually these should likely all become ForceProviders.
563  * Within this function the intention is to have algorithms that do
564  * global communication at the end, so global barriers within the MD loop
565  * are as close together as possible.
566  *
567  * \param[in]     fplog            The log file
568  * \param[in]     cr               The communication record
569  * \param[in]     inputrec         The input record
570  * \param[in]     awh              The Awh module (nullptr if none in use).
571  * \param[in]     enforcedRotation Enforced rotation module.
572  * \param[in]     imdSession       The IMD session
573  * \param[in]     pull_work        The pull work structure.
574  * \param[in]     step             The current MD step
575  * \param[in]     t                The current time
576  * \param[in,out] wcycle           Wallcycle accounting struct
577  * \param[in,out] forceProviders   Pointer to a list of force providers
578  * \param[in]     box              The unit cell
579  * \param[in]     x                The coordinates
580  * \param[in]     mdatoms          Per atom properties
581  * \param[in]     lambda           Array of free-energy lambda values
582  * \param[in]     stepWork         Step schedule flags
583  * \param[in,out] forceWithVirial  Force and virial buffers
584  * \param[in,out] enerd            Energy buffer
585  * \param[in,out] ed               Essential dynamics pointer
586  * \param[in]     didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
587  *
588  * \todo Remove didNeighborSearch, which is used incorrectly.
589  * \todo Convert all other algorithms called here to ForceProviders.
590  */
591 static void computeSpecialForces(FILE*                          fplog,
592                                  const t_commrec*               cr,
593                                  const t_inputrec*              inputrec,
594                                  gmx::Awh*                      awh,
595                                  gmx_enfrot*                    enforcedRotation,
596                                  gmx::ImdSession*               imdSession,
597                                  pull_t*                        pull_work,
598                                  int64_t                        step,
599                                  double                         t,
600                                  gmx_wallcycle_t                wcycle,
601                                  gmx::ForceProviders*           forceProviders,
602                                  const matrix                   box,
603                                  gmx::ArrayRef<const gmx::RVec> x,
604                                  const t_mdatoms*               mdatoms,
605                                  gmx::ArrayRef<const real>      lambda,
606                                  const StepWorkload&            stepWork,
607                                  gmx::ForceWithVirial*          forceWithVirial,
608                                  gmx_enerdata_t*                enerd,
609                                  gmx_edsam*                     ed,
610                                  bool                           didNeighborSearch)
611 {
612     /* NOTE: Currently all ForceProviders only provide forces.
613      *       When they also provide energies, remove this conditional.
614      */
615     if (stepWork.computeForces)
616     {
617         gmx::ForceProviderInput  forceProviderInput(x, *mdatoms, t, box, *cr);
618         gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
619
620         /* Collect forces from modules */
621         forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
622     }
623
624     if (inputrec->bPull && pull_have_potential(pull_work))
625     {
626         pull_potential_wrapper(cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work,
627                                lambda.data(), t, wcycle);
628     }
629     if (awh)
630     {
631         const bool          needForeignEnergyDifferences = awh->needForeignEnergyDifferences(step);
632         std::vector<double> foreignLambdaDeltaH, foreignLambdaDhDl;
633         if (needForeignEnergyDifferences)
634         {
635             enerd->foreignLambdaTerms.finalizePotentialContributions(enerd->dvdl_lin, lambda,
636                                                                      *inputrec->fepvals);
637             std::tie(foreignLambdaDeltaH, foreignLambdaDhDl) = enerd->foreignLambdaTerms.getTerms(cr);
638         }
639
640         enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(
641                 inputrec->pbcType, mdatoms->massT, foreignLambdaDeltaH, foreignLambdaDhDl, box,
642                 forceWithVirial, t, step, wcycle, fplog);
643     }
644
645     rvec* f = as_rvec_array(forceWithVirial->force_.data());
646
647     /* Add the forces from enforced rotation potentials (if any) */
648     if (inputrec->bRot)
649     {
650         wallcycle_start(wcycle, ewcROTadd);
651         enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
652         wallcycle_stop(wcycle, ewcROTadd);
653     }
654
655     if (ed)
656     {
657         /* Note that since init_edsam() is called after the initialization
658          * of forcerec, edsam doesn't request the noVirSum force buffer.
659          * Thus if no other algorithm (e.g. PME) requires it, the forces
660          * here will contribute to the virial.
661          */
662         do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
663     }
664
665     /* Add forces from interactive molecular dynamics (IMD), if any */
666     if (inputrec->bIMD && stepWork.computeForces)
667     {
668         imdSession->applyForces(f);
669     }
670 }
671
672 /*! \brief Launch the prepare_step and spread stages of PME GPU.
673  *
674  * \param[in]  pmedata              The PME structure
675  * \param[in]  box                  The box matrix
676  * \param[in]  stepWork             Step schedule flags
677  * \param[in]  xReadyOnDevice       Event synchronizer indicating that the coordinates are ready in the device memory.
678  * \param[in]  lambdaQ              The Coulomb lambda of the current state.
679  * \param[in]  wcycle               The wallcycle structure
680  */
681 static inline void launchPmeGpuSpread(gmx_pme_t*            pmedata,
682                                       const matrix          box,
683                                       const StepWorkload&   stepWork,
684                                       GpuEventSynchronizer* xReadyOnDevice,
685                                       const real            lambdaQ,
686                                       gmx_wallcycle_t       wcycle)
687 {
688     pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
689     pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle, lambdaQ);
690 }
691
692 /*! \brief Launch the FFT and gather stages of PME GPU
693  *
694  * This function only implements setting the output forces (no accumulation).
695  *
696  * \param[in]  pmedata        The PME structure
697  * \param[in]  lambdaQ        The Coulomb lambda of the current system state.
698  * \param[in]  wcycle         The wallcycle structure
699  * \param[in]  stepWork       Step schedule flags
700  */
701 static void launchPmeGpuFftAndGather(gmx_pme_t*               pmedata,
702                                      const real               lambdaQ,
703                                      gmx_wallcycle_t          wcycle,
704                                      const gmx::StepWorkload& stepWork)
705 {
706     pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
707     pme_gpu_launch_gather(pmedata, wcycle, lambdaQ);
708 }
709
710 /*! \brief
711  *  Polling wait for either of the PME or nonbonded GPU tasks.
712  *
713  * Instead of a static order in waiting for GPU tasks, this function
714  * polls checking which of the two tasks completes first, and does the
715  * associated force buffer reduction overlapped with the other task.
716  * By doing that, unlike static scheduling order, it can always overlap
717  * one of the reductions, regardless of the GPU task completion order.
718  *
719  * \param[in]     nbv              Nonbonded verlet structure
720  * \param[in,out] pmedata          PME module data
721  * \param[in,out] forceOutputsNonbonded  Force outputs for the non-bonded forces and shift forces
722  * \param[in,out] forceOutputsPme  Force outputs for the PME forces and virial
723  * \param[in,out] enerd            Energy data structure results are reduced into
724  * \param[in]     lambdaQ          The Coulomb lambda of the current system state.
725  * \param[in]     stepWork         Step schedule flags
726  * \param[in]     wcycle           The wallcycle structure
727  */
728 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
729                                         gmx_pme_t*          pmedata,
730                                         gmx::ForceOutputs*  forceOutputsNonbonded,
731                                         gmx::ForceOutputs*  forceOutputsPme,
732                                         gmx_enerdata_t*     enerd,
733                                         const real          lambdaQ,
734                                         const StepWorkload& stepWork,
735                                         gmx_wallcycle_t     wcycle)
736 {
737     bool isPmeGpuDone = false;
738     bool isNbGpuDone  = false;
739
740     gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
741
742     while (!isPmeGpuDone || !isNbGpuDone)
743     {
744         if (!isPmeGpuDone)
745         {
746             GpuTaskCompletion completionType =
747                     (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
748             isPmeGpuDone = pme_gpu_try_finish_task(pmedata, stepWork, wcycle,
749                                                    &forceOutputsPme->forceWithVirial(), enerd,
750                                                    lambdaQ, completionType);
751         }
752
753         if (!isNbGpuDone)
754         {
755             auto&             forceBuffersNonbonded = forceOutputsNonbonded->forceWithShiftForces();
756             GpuTaskCompletion completionType =
757                     (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
758             isNbGpuDone = Nbnxm::gpu_try_finish_task(
759                     nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
760                     enerd->grpp.ener[egCOULSR].data(), forceBuffersNonbonded.shiftForces(),
761                     completionType, wcycle);
762
763             if (isNbGpuDone)
764             {
765                 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceBuffersNonbonded.force());
766             }
767         }
768     }
769 }
770
771 /*! \brief Set up the different force buffers; also does clearing.
772  *
773  * \param[in] forceHelperBuffers  Helper force buffers
774  * \param[in] force     force array
775  * \param[in] stepWork  Step schedule flags
776  * \param[out] wcycle   wallcycle recording structure
777  *
778  * \returns             Cleared force output structure
779  */
780 static ForceOutputs setupForceOutputs(ForceHelperBuffers*                 forceHelperBuffers,
781                                       gmx::ArrayRefWithPadding<gmx::RVec> force,
782                                       const StepWorkload&                 stepWork,
783                                       gmx_wallcycle_t                     wcycle)
784 {
785     wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
786
787     /* NOTE: We assume fr->shiftForces is all zeros here */
788     gmx::ForceWithShiftForces forceWithShiftForces(force, stepWork.computeVirial,
789                                                    forceHelperBuffers->shiftForces());
790
791     if (stepWork.computeForces)
792     {
793         /* Clear the short- and long-range forces */
794         clearRVecs(forceWithShiftForces.force(), true);
795
796         /* Clear the shift forces */
797         clearRVecs(forceWithShiftForces.shiftForces(), false);
798     }
799
800     /* If we need to compute the virial, we might need a separate
801      * force buffer for algorithms for which the virial is calculated
802      * directly, such as PME. Otherwise, forceWithVirial uses the
803      * the same force (f in legacy calls) buffer as other algorithms.
804      */
805     const bool useSeparateForceWithVirialBuffer =
806             (stepWork.computeForces
807              && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
808     /* forceWithVirial uses the local atom range only */
809     gmx::ForceWithVirial forceWithVirial(
810             useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
811                                              : force.unpaddedArrayRef(),
812             stepWork.computeVirial);
813
814     if (useSeparateForceWithVirialBuffer)
815     {
816         /* TODO: update comment
817          * We only compute forces on local atoms. Note that vsites can
818          * spread to non-local atoms, but that part of the buffer is
819          * cleared separately in the vsite spreading code.
820          */
821         clearRVecs(forceWithVirial.force_, true);
822     }
823
824     wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
825
826     return ForceOutputs(forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(),
827                         forceWithVirial);
828 }
829
830
831 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
832  */
833 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec&         inputrec,
834                                                           const t_forcerec&         fr,
835                                                           const pull_t*             pull_work,
836                                                           const gmx_edsam*          ed,
837                                                           const t_mdatoms&          mdatoms,
838                                                           const SimulationWorkload& simulationWork,
839                                                           const StepWorkload&       stepWork)
840 {
841     DomainLifetimeWorkload domainWork;
842     // Note that haveSpecialForces is constant over the whole run
843     domainWork.haveSpecialForces =
844             haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
845     domainWork.haveCpuListedForceWork = false;
846     domainWork.haveCpuBondedWork      = false;
847     for (const auto& listedForces : fr.listedForces)
848     {
849         if (listedForces.haveCpuListedForces(*fr.fcdata))
850         {
851             domainWork.haveCpuListedForceWork = true;
852         }
853         if (listedForces.haveCpuBondeds())
854         {
855             domainWork.haveCpuBondedWork = true;
856         }
857     }
858     domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
859     // Note that haveFreeEnergyWork is constant over the whole run
860     domainWork.haveFreeEnergyWork = (fr.efep != efepNO && mdatoms.nPerturbed != 0);
861     // We assume we have local force work if there are CPU
862     // force tasks including PME or nonbondeds.
863     domainWork.haveCpuLocalForceWork =
864             domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
865             || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
866             || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
867
868     return domainWork;
869 }
870
871 /*! \brief Set up force flag stuct from the force bitmask.
872  *
873  * \param[in]      legacyFlags          Force bitmask flags used to construct the new flags
874  * \param[in]      mtsLevels            The multiple time-stepping levels, either empty or 2 levels
875  * \param[in]      step                 The current MD step
876  * \param[in]      simulationWork       Simulation workload description.
877  * \param[in]      rankHasPmeDuty       If this rank computes PME.
878  *
879  * \returns New Stepworkload description.
880  */
881 static StepWorkload setupStepWorkload(const int                     legacyFlags,
882                                       ArrayRef<const gmx::MtsLevel> mtsLevels,
883                                       const int64_t                 step,
884                                       const SimulationWorkload&     simulationWork,
885                                       const bool                    rankHasPmeDuty)
886 {
887     GMX_ASSERT(mtsLevels.empty() || mtsLevels.size() == 2, "Expect 0 or 2 MTS levels");
888     const bool computeSlowForces = (mtsLevels.empty() || step % mtsLevels[1].stepFactor == 0);
889
890     StepWorkload flags;
891     flags.stateChanged        = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
892     flags.haveDynamicBox      = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
893     flags.doNeighborSearch    = ((legacyFlags & GMX_FORCE_NS) != 0);
894     flags.computeSlowForces   = computeSlowForces;
895     flags.computeVirial       = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
896     flags.computeEnergy       = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
897     flags.computeForces       = ((legacyFlags & GMX_FORCE_FORCES) != 0);
898     flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
899     flags.computeNonbondedForces =
900             ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded
901             && !(simulationWork.computeNonbondedAtMtsLevel1 && !computeSlowForces);
902     flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
903
904     if (simulationWork.useGpuBufferOps)
905     {
906         GMX_ASSERT(simulationWork.useGpuNonbonded,
907                    "Can only offload buffer ops if nonbonded computation is also offloaded");
908     }
909     flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
910     // on virial steps the CPU reduction path is taken
911     flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
912     flags.useGpuPmeFReduction = flags.computeSlowForces && flags.useGpuFBufferOps && simulationWork.useGpuPme
913                                 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication);
914     flags.useGpuXHalo = simulationWork.useGpuHaloExchange;
915     flags.useGpuFHalo = simulationWork.useGpuHaloExchange && flags.useGpuFBufferOps;
916
917     return flags;
918 }
919
920
921 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
922  *
923  * TODO: eliminate \p useGpuPmeOnThisRank when this is
924  * incorporated in DomainLifetimeWorkload.
925  */
926 static void launchGpuEndOfStepTasks(nonbonded_verlet_t*               nbv,
927                                     gmx::GpuBonded*                   gpuBonded,
928                                     gmx_pme_t*                        pmedata,
929                                     gmx_enerdata_t*                   enerd,
930                                     const gmx::MdrunScheduleWorkload& runScheduleWork,
931                                     bool                              useGpuPmeOnThisRank,
932                                     int64_t                           step,
933                                     gmx_wallcycle_t                   wcycle)
934 {
935     if (runScheduleWork.simulationWork.useGpuNonbonded && runScheduleWork.stepWork.computeNonbondedForces)
936     {
937         /* Launch pruning before buffer clearing because the API overhead of the
938          * clear kernel launches can leave the GPU idle while it could be running
939          * the prune kernel.
940          */
941         if (nbv->isDynamicPruningStepGpu(step))
942         {
943             nbv->dispatchPruneKernelGpu(step);
944         }
945
946         /* now clear the GPU outputs while we finish the step on the CPU */
947         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
948         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
949         Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
950         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
951         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
952     }
953
954     if (useGpuPmeOnThisRank)
955     {
956         pme_gpu_reinit_computation(pmedata, wcycle);
957     }
958
959     if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
960     {
961         // in principle this should be included in the DD balancing region,
962         // but generally it is infrequent so we'll omit it for the sake of
963         // simpler code
964         gpuBonded->waitAccumulateEnergyTerms(enerd);
965
966         gpuBonded->clearEnergies();
967     }
968 }
969
970 //! \brief Data structure to hold dipole-related data and staging arrays
971 struct DipoleData
972 {
973     //! Dipole staging for fast summing over MPI
974     gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
975     //! Dipole staging for states A and B (index 0 and 1 resp.)
976     gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
977 };
978
979
980 static void reduceAndUpdateMuTot(DipoleData*                   dipoleData,
981                                  const t_commrec*              cr,
982                                  const bool                    haveFreeEnergy,
983                                  gmx::ArrayRef<const real>     lambda,
984                                  rvec                          muTotal,
985                                  const DDBalanceRegionHandler& ddBalanceRegionHandler)
986 {
987     if (PAR(cr))
988     {
989         gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
990         ddBalanceRegionHandler.reopenRegionCpu();
991     }
992     for (int i = 0; i < 2; i++)
993     {
994         for (int j = 0; j < DIM; j++)
995         {
996             dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
997         }
998     }
999
1000     if (!haveFreeEnergy)
1001     {
1002         copy_rvec(dipoleData->muStateAB[0], muTotal);
1003     }
1004     else
1005     {
1006         for (int j = 0; j < DIM; j++)
1007         {
1008             muTotal[j] = (1.0 - lambda[efptCOUL]) * dipoleData->muStateAB[0][j]
1009                          + lambda[efptCOUL] * dipoleData->muStateAB[1][j];
1010         }
1011     }
1012 }
1013
1014 /*! \brief Combines MTS level0 and level1 force buffes into a full and MTS-combined force buffer.
1015  *
1016  * \param[in]     numAtoms        The number of atoms to combine forces for
1017  * \param[in,out] forceMtsLevel0  Input: F_level0, output: F_level0 + F_level1
1018  * \param[in,out] forceMts        Input: F_level1, output: F_level0 + mtsFactor * F_level1
1019  * \param[in]     mtsFactor       The factor between the level0 and level1 time step
1020  */
1021 static void combineMtsForces(const int      numAtoms,
1022                              ArrayRef<RVec> forceMtsLevel0,
1023                              ArrayRef<RVec> forceMts,
1024                              const real     mtsFactor)
1025 {
1026     const int gmx_unused numThreads = gmx_omp_nthreads_get(emntDefault);
1027 #pragma omp parallel for num_threads(numThreads) schedule(static)
1028     for (int i = 0; i < numAtoms; i++)
1029     {
1030         const RVec forceMtsLevel0Tmp = forceMtsLevel0[i];
1031         forceMtsLevel0[i] += forceMts[i];
1032         forceMts[i] = forceMtsLevel0Tmp + mtsFactor * forceMts[i];
1033     }
1034 }
1035
1036 /*! \brief Setup for the local and non-local GPU force reductions:
1037  * reinitialization plus the registration of forces and dependencies.
1038  *
1039  * \param [in] runScheduleWork               Schedule workload flag structure
1040  * \param [in] cr                            Communication record object
1041  * \param [in] fr                            Force record object
1042  */
1043 static void setupGpuForceReductions(gmx::MdrunScheduleWorkload* runScheduleWork,
1044                                     const t_commrec*            cr,
1045                                     t_forcerec*                 fr)
1046 {
1047
1048     nonbonded_verlet_t*          nbv      = fr->nbv.get();
1049     gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1050
1051     // (re-)initialize local GPU force reduction
1052     const bool accumulate =
1053             runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr);
1054     const int atomStart = 0;
1055     fr->gpuForceReduction[gmx::AtomLocality::Local]->reinit(
1056             stateGpu->getForces(), nbv->getNumAtoms(AtomLocality::Local), nbv->getGridIndices(),
1057             atomStart, accumulate, stateGpu->fReducedOnDevice());
1058
1059     // register forces and add dependencies
1060     fr->gpuForceReduction[gmx::AtomLocality::Local]->registerNbnxmForce(nbv->getGpuForces());
1061
1062     if (runScheduleWork->simulationWork.useGpuPme
1063         && (thisRankHasDuty(cr, DUTY_PME) || runScheduleWork->simulationWork.useGpuPmePpCommunication))
1064     {
1065         void* forcePtr = thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1066                                                        : // PME force buffer on same GPU
1067                                  fr->pmePpCommGpu->getGpuForceStagingPtr(); // buffer received from other GPU
1068         fr->gpuForceReduction[gmx::AtomLocality::Local]->registerRvecForce(forcePtr);
1069
1070         GpuEventSynchronizer* const pmeSynchronizer =
1071                 (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1072                                                : // PME force buffer on same GPU
1073                          fr->pmePpCommGpu->getForcesReadySynchronizer()); // buffer received from other GPU
1074         fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(pmeSynchronizer);
1075     }
1076
1077     if ((runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr))
1078         && !runScheduleWork->simulationWork.useGpuHaloExchange)
1079     {
1080         fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1081                 stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::Local, true));
1082     }
1083
1084     if (runScheduleWork->simulationWork.useGpuHaloExchange)
1085     {
1086         fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1087                 cr->dd->gpuHaloExchange[0][0]->getForcesReadyOnDeviceEvent());
1088     }
1089
1090     if (havePPDomainDecomposition(cr))
1091     {
1092         // (re-)initialize non-local GPU force reduction
1093         const bool accumulate = runScheduleWork->domainWork.haveCpuBondedWork
1094                                 || runScheduleWork->domainWork.haveFreeEnergyWork;
1095         const int atomStart = dd_numHomeAtoms(*cr->dd);
1096         fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->reinit(
1097                 stateGpu->getForces(), nbv->getNumAtoms(AtomLocality::NonLocal),
1098                 nbv->getGridIndices(), atomStart, accumulate);
1099
1100         // register forces and add dependencies
1101         fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->registerNbnxmForce(nbv->getGpuForces());
1102         if (runScheduleWork->domainWork.haveCpuBondedWork || runScheduleWork->domainWork.haveFreeEnergyWork)
1103         {
1104             fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->addDependency(
1105                     stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::NonLocal, true));
1106         }
1107     }
1108 }
1109
1110
1111 void do_force(FILE*                               fplog,
1112               const t_commrec*                    cr,
1113               const gmx_multisim_t*               ms,
1114               const t_inputrec*                   inputrec,
1115               gmx::Awh*                           awh,
1116               gmx_enfrot*                         enforcedRotation,
1117               gmx::ImdSession*                    imdSession,
1118               pull_t*                             pull_work,
1119               int64_t                             step,
1120               t_nrnb*                             nrnb,
1121               gmx_wallcycle_t                     wcycle,
1122               const gmx_localtop_t*               top,
1123               const matrix                        box,
1124               gmx::ArrayRefWithPadding<gmx::RVec> x,
1125               history_t*                          hist,
1126               gmx::ForceBuffersView*              forceView,
1127               tensor                              vir_force,
1128               const t_mdatoms*                    mdatoms,
1129               gmx_enerdata_t*                     enerd,
1130               gmx::ArrayRef<const real>           lambda,
1131               t_forcerec*                         fr,
1132               gmx::MdrunScheduleWorkload*         runScheduleWork,
1133               gmx::VirtualSitesHandler*           vsite,
1134               rvec                                muTotal,
1135               double                              t,
1136               gmx_edsam*                          ed,
1137               int                                 legacyFlags,
1138               const DDBalanceRegionHandler&       ddBalanceRegionHandler)
1139 {
1140     auto force = forceView->forceWithPadding();
1141     GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1142                "The size of the force buffer should be at least the number of atoms to compute "
1143                "forces for");
1144
1145     nonbonded_verlet_t*  nbv = fr->nbv.get();
1146     interaction_const_t* ic  = fr->ic;
1147
1148     gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1149
1150     const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1151
1152     runScheduleWork->stepWork    = setupStepWorkload(legacyFlags, inputrec->mtsLevels, step,
1153                                                   simulationWork, thisRankHasDuty(cr, DUTY_PME));
1154     const StepWorkload& stepWork = runScheduleWork->stepWork;
1155
1156     const bool useGpuPmeOnThisRank =
1157             simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces;
1158
1159     /* At a search step we need to start the first balancing region
1160      * somewhere early inside the step after communication during domain
1161      * decomposition (and not during the previous step as usual).
1162      */
1163     if (stepWork.doNeighborSearch)
1164     {
1165         ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1166     }
1167
1168     clear_mat(vir_force);
1169
1170     if (fr->pbcType != PbcType::No)
1171     {
1172         /* Compute shift vectors every step,
1173          * because of pressure coupling or box deformation!
1174          */
1175         if (stepWork.haveDynamicBox && stepWork.stateChanged)
1176         {
1177             calc_shifts(box, fr->shift_vec);
1178         }
1179
1180         const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1181         const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1182         if (calcCGCM)
1183         {
1184             put_atoms_in_box_omp(fr->pbcType, box, x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1185                                  gmx_omp_nthreads_get(emntDefault));
1186             inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1187         }
1188     }
1189
1190     nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1191
1192     const bool pmeSendCoordinatesFromGpu =
1193             GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1194     const bool reinitGpuPmePpComms =
1195             GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1196
1197     const auto localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1198                                              ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1199                                                        AtomLocality::Local, simulationWork, stepWork)
1200                                              : nullptr;
1201
1202     // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1203     // Otherwise the send will occur after H2D coordinate transfer.
1204     if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu && stepWork.computeSlowForces)
1205     {
1206         /* Send particle coordinates to the pme nodes */
1207         if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1208         {
1209             GMX_RELEASE_ASSERT(false,
1210                                "GPU update and separate PME ranks are only supported with GPU "
1211                                "direct communication!");
1212             // TODO: when this code-path becomes supported add:
1213             // stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1214         }
1215
1216         gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1217                                  lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1218                                  step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1219                                  pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1220     }
1221
1222     // Coordinates on the device are needed if PME or BufferOps are offloaded.
1223     // The local coordinates can be copied right away.
1224     // NOTE: Consider moving this copy to right after they are updated and constrained,
1225     //       if the later is not offloaded.
1226     if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1227     {
1228         if (stepWork.doNeighborSearch)
1229         {
1230             // TODO refactor this to do_md, after partitioning.
1231             stateGpu->reinit(mdatoms->homenr,
1232                              cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1233             if (useGpuPmeOnThisRank)
1234             {
1235                 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1236                 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1237             }
1238         }
1239         // We need to copy coordinates when:
1240         // 1. Update is not offloaded
1241         // 2. The buffers were reinitialized on search step
1242         if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1243         {
1244             GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1245             stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1246         }
1247     }
1248
1249     // Copy coordinate from the GPU if update is on the GPU and there
1250     // are forces to be computed on the CPU, or for the computation of
1251     // virial, or if host-side data will be transferred from this task
1252     // to a remote task for halo exchange or PME-PP communication. At
1253     // search steps the current coordinates are already on the host,
1254     // hence copy is not needed.
1255     const bool haveHostPmePpComms =
1256             !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1257
1258     GMX_ASSERT(simulationWork.useGpuHaloExchange
1259                        == ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange[0].empty())),
1260                "The GPU halo exchange is active, but it has not been constructed.");
1261     const bool haveHostHaloExchangeComms =
1262             havePPDomainDecomposition(cr) && !simulationWork.useGpuHaloExchange;
1263
1264     bool gmx_used_in_debug haveCopiedXFromGpu = false;
1265     if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1266         && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1267             || haveHostPmePpComms || haveHostHaloExchangeComms))
1268     {
1269         GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1270         stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1271         haveCopiedXFromGpu = true;
1272     }
1273
1274     // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1275     // Otherwise the send will occur before the H2D coordinate transfer.
1276     if (!thisRankHasDuty(cr, DUTY_PME) && pmeSendCoordinatesFromGpu)
1277     {
1278         /* Send particle coordinates to the pme nodes */
1279         gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1280                                  lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1281                                  step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1282                                  pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1283     }
1284
1285     if (useGpuPmeOnThisRank)
1286     {
1287         launchPmeGpuSpread(fr->pmedata, box, stepWork, localXReadyOnDevice, lambda[efptCOUL], wcycle);
1288     }
1289
1290     const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1291
1292     /* do gridding for pair search */
1293     if (stepWork.doNeighborSearch)
1294     {
1295         if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1296         {
1297             fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1298         }
1299
1300         // TODO
1301         // - vzero is constant, do we need to pass it?
1302         // - box_diag should be passed directly to nbnxn_put_on_grid
1303         //
1304         rvec vzero;
1305         clear_rvec(vzero);
1306
1307         rvec box_diag;
1308         box_diag[XX] = box[XX][XX];
1309         box_diag[YY] = box[YY][YY];
1310         box_diag[ZZ] = box[ZZ][ZZ];
1311
1312         wallcycle_start(wcycle, ewcNS);
1313         if (!DOMAINDECOMP(cr))
1314         {
1315             wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1316             nbnxn_put_on_grid(nbv, box, 0, vzero, box_diag, nullptr, { 0, mdatoms->homenr }, -1,
1317                               fr->cginfo, x.unpaddedArrayRef(), 0, nullptr);
1318             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1319         }
1320         else
1321         {
1322             wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1323             nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1324             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1325         }
1326
1327         nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1328                                gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr), fr->cginfo);
1329
1330         wallcycle_stop(wcycle, ewcNS);
1331
1332         /* initialize the GPU nbnxm atom data and bonded data structures */
1333         if (simulationWork.useGpuNonbonded)
1334         {
1335             wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1336
1337             wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1338             Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1339             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1340
1341             if (fr->gpuBonded)
1342             {
1343                 /* Now we put all atoms on the grid, we can assign bonded
1344                  * interactions to the GPU, where the grid order is
1345                  * needed. Also the xq, f and fshift device buffers have
1346                  * been reallocated if needed, so the bonded code can
1347                  * learn about them. */
1348                 // TODO the xq, f, and fshift buffers are now shared
1349                 // resources, so they should be maintained by a
1350                 // higher-level object than the nb module.
1351                 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(
1352                         nbv->getGridIndices(), top->idef, Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1353                         Nbnxm::gpu_get_f(nbv->gpu_nbv), Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1354             }
1355             wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1356         }
1357
1358         // Need to run after the GPU-offload bonded interaction lists
1359         // are set up to be able to determine whether there is bonded work.
1360         runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1361                 *inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1362
1363         wallcycle_start_nocount(wcycle, ewcNS);
1364         wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1365         /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1366         nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1367
1368         nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1369
1370         wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1371         wallcycle_stop(wcycle, ewcNS);
1372
1373         if (stepWork.useGpuXBufferOps)
1374         {
1375             nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1376         }
1377
1378         if (simulationWork.useGpuBufferOps)
1379         {
1380             setupGpuForceReductions(runScheduleWork, cr, fr);
1381         }
1382     }
1383     else if (!EI_TPI(inputrec->eI) && stepWork.computeNonbondedForces)
1384     {
1385         if (stepWork.useGpuXBufferOps)
1386         {
1387             GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1388             nbv->convertCoordinatesGpu(AtomLocality::Local, false, stateGpu->getCoordinates(),
1389                                        localXReadyOnDevice);
1390         }
1391         else
1392         {
1393             if (simulationWork.useGpuUpdate)
1394             {
1395                 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1396                 GMX_ASSERT(haveCopiedXFromGpu,
1397                            "a wait should only be triggered if copy has been scheduled");
1398                 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1399             }
1400             nbv->convertCoordinates(AtomLocality::Local, false, x.unpaddedArrayRef());
1401         }
1402     }
1403
1404     if (simulationWork.useGpuNonbonded && (stepWork.computeNonbondedForces || domainWork.haveGpuBondedWork))
1405     {
1406         ddBalanceRegionHandler.openBeforeForceComputationGpu();
1407
1408         wallcycle_start(wcycle, ewcLAUNCH_GPU);
1409
1410         wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1411         Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1412         if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1413         {
1414             Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1415         }
1416         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1417         // with X buffer ops offloaded to the GPU on all but the search steps
1418
1419         // bonded work not split into separate local and non-local, so with DD
1420         // we can only launch the kernel after non-local coordinates have been received.
1421         if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1422         {
1423             wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1424             fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1425             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1426         }
1427
1428         /* launch local nonbonded work on GPU */
1429         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1430         do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1431         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1432         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1433     }
1434
1435     if (useGpuPmeOnThisRank)
1436     {
1437         // In PME GPU and mixed mode we launch FFT / gather after the
1438         // X copy/transform to allow overlap as well as after the GPU NB
1439         // launch to avoid FFT launch overhead hijacking the CPU and delaying
1440         // the nonbonded kernel.
1441         launchPmeGpuFftAndGather(fr->pmedata, lambda[efptCOUL], wcycle, stepWork);
1442     }
1443
1444     /* Communicate coordinates and sum dipole if necessary +
1445        do non-local pair search */
1446     if (havePPDomainDecomposition(cr))
1447     {
1448         if (stepWork.doNeighborSearch)
1449         {
1450             // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1451             wallcycle_start_nocount(wcycle, ewcNS);
1452             wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1453             /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1454             nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1455
1456             nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1457             wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1458             wallcycle_stop(wcycle, ewcNS);
1459             // TODO refactor this GPU halo exchange re-initialisation
1460             // to location in do_md where GPU halo exchange is
1461             // constructed at partitioning, after above stateGpu
1462             // re-initialization has similarly been refactored
1463             if (simulationWork.useGpuHaloExchange)
1464             {
1465                 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1466             }
1467         }
1468         else
1469         {
1470             if (stepWork.useGpuXHalo)
1471             {
1472                 // The following must be called after local setCoordinates (which records an event
1473                 // when the coordinate data has been copied to the device).
1474                 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1475
1476                 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1477                 {
1478                     // non-local part of coordinate buffer must be copied back to host for CPU work
1479                     stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1480                 }
1481             }
1482             else
1483             {
1484                 // Note: GPU update + DD without direct communication is not supported,
1485                 // a waitCoordinatesReadyOnHost() should be issued if it will be.
1486                 GMX_ASSERT(!simulationWork.useGpuUpdate,
1487                            "GPU update is not supported with CPU halo exchange");
1488                 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1489             }
1490
1491             if (stepWork.useGpuXBufferOps)
1492             {
1493                 if (!useGpuPmeOnThisRank && !stepWork.useGpuXHalo)
1494                 {
1495                     stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1496                 }
1497                 nbv->convertCoordinatesGpu(AtomLocality::NonLocal, false, stateGpu->getCoordinates(),
1498                                            stateGpu->getCoordinatesReadyOnDeviceEvent(
1499                                                    AtomLocality::NonLocal, simulationWork, stepWork));
1500             }
1501             else
1502             {
1503                 nbv->convertCoordinates(AtomLocality::NonLocal, false, x.unpaddedArrayRef());
1504             }
1505         }
1506
1507         if (simulationWork.useGpuNonbonded)
1508         {
1509             wallcycle_start(wcycle, ewcLAUNCH_GPU);
1510
1511             if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1512             {
1513                 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1514                 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1515                 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1516             }
1517
1518             if (domainWork.haveGpuBondedWork)
1519             {
1520                 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1521                 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1522                 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1523             }
1524
1525             /* launch non-local nonbonded tasks on GPU */
1526             wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1527             do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1528                          nrnb, wcycle);
1529             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1530
1531             wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1532         }
1533     }
1534
1535     if (simulationWork.useGpuNonbonded && stepWork.computeNonbondedForces)
1536     {
1537         /* launch D2H copy-back F */
1538         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1539         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1540
1541         if (havePPDomainDecomposition(cr))
1542         {
1543             Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1544         }
1545         Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1546         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1547
1548         if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1549         {
1550             fr->gpuBonded->launchEnergyTransfer();
1551         }
1552         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1553     }
1554
1555     gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1556     if (fr->wholeMoleculeTransform)
1557     {
1558         xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1559     }
1560
1561     DipoleData dipoleData;
1562
1563     if (simulationWork.computeMuTot)
1564     {
1565         const int start = 0;
1566
1567         /* Calculate total (local) dipole moment in a temporary common array.
1568          * This makes it possible to sum them over nodes faster.
1569          */
1570         gmx::ArrayRef<const gmx::RVec> xRef =
1571                 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1572         calc_mu(start, mdatoms->homenr, xRef, mdatoms->chargeA, mdatoms->chargeB,
1573                 mdatoms->nChargePerturbed, dipoleData.muStaging[0], dipoleData.muStaging[1]);
1574
1575         reduceAndUpdateMuTot(&dipoleData, cr, (fr->efep != efepNO), lambda, muTotal, ddBalanceRegionHandler);
1576     }
1577
1578     /* Reset energies */
1579     reset_enerdata(enerd);
1580
1581     if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1582     {
1583         wallcycle_start(wcycle, ewcPPDURINGPME);
1584         dd_force_flop_start(cr->dd, nrnb);
1585     }
1586
1587     // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1588     // this wait ensures that the D2H transfer is complete.
1589     if ((simulationWork.useGpuUpdate)
1590         && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1591     {
1592         stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1593     }
1594
1595     if (inputrec->bRot)
1596     {
1597         wallcycle_start(wcycle, ewcROT);
1598         do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step,
1599                     stepWork.doNeighborSearch);
1600         wallcycle_stop(wcycle, ewcROT);
1601     }
1602
1603     /* Start the force cycle counter.
1604      * Note that a different counter is used for dynamic load balancing.
1605      */
1606     wallcycle_start(wcycle, ewcFORCE);
1607
1608     /* Set up and clear force outputs:
1609      * forceOutMtsLevel0:  everything except what is in the other two outputs
1610      * forceOutMtsLevel1:  PME-mesh and listed-forces group 1
1611      * forceOutNonbonded: non-bonded forces
1612      * Without multiple time stepping all point to the same object.
1613      * With multiple time-stepping the use is different for MTS fast (level0 only) and slow steps.
1614      */
1615     ForceOutputs forceOutMtsLevel0 =
1616             setupForceOutputs(&fr->forceHelperBuffers[0], force, stepWork, wcycle);
1617
1618     // Force output for MTS combined forces, only set at level1 MTS steps
1619     std::optional<ForceOutputs> forceOutMts =
1620             (fr->useMts && stepWork.computeSlowForces)
1621                     ? std::optional(setupForceOutputs(&fr->forceHelperBuffers[1],
1622                                                       forceView->forceMtsCombinedWithPadding(),
1623                                                       stepWork, wcycle))
1624                     : std::nullopt;
1625
1626     ForceOutputs* forceOutMtsLevel1 =
1627             fr->useMts ? (stepWork.computeSlowForces ? &forceOutMts.value() : nullptr) : &forceOutMtsLevel0;
1628
1629     const bool nonbondedAtMtsLevel1 = runScheduleWork->simulationWork.computeNonbondedAtMtsLevel1;
1630
1631     ForceOutputs* forceOutNonbonded = nonbondedAtMtsLevel1 ? forceOutMtsLevel1 : &forceOutMtsLevel0;
1632
1633     if (inputrec->bPull && pull_have_constraint(pull_work))
1634     {
1635         clear_pull_forces(pull_work);
1636     }
1637
1638     /* We calculate the non-bonded forces, when done on the CPU, here.
1639      * We do this before calling do_force_lowlevel, because in that
1640      * function, the listed forces are calculated before PME, which
1641      * does communication.  With this order, non-bonded and listed
1642      * force calculation imbalance can be balanced out by the domain
1643      * decomposition load balancing.
1644      */
1645
1646     const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1647
1648     if (!useOrEmulateGpuNb)
1649     {
1650         do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1651     }
1652
1653     if (fr->efep != efepNO && stepWork.computeNonbondedForces)
1654     {
1655         /* Calculate the local and non-local free energy interactions here.
1656          * Happens here on the CPU both with and without GPU.
1657          */
1658         nbv->dispatchFreeEnergyKernel(InteractionLocality::Local, fr,
1659                                       as_rvec_array(x.unpaddedArrayRef().data()),
1660                                       &forceOutNonbonded->forceWithShiftForces(), *mdatoms,
1661                                       inputrec->fepvals, lambda, enerd, stepWork, nrnb);
1662
1663         if (havePPDomainDecomposition(cr))
1664         {
1665             nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal, fr,
1666                                           as_rvec_array(x.unpaddedArrayRef().data()),
1667                                           &forceOutNonbonded->forceWithShiftForces(), *mdatoms,
1668                                           inputrec->fepvals, lambda, enerd, stepWork, nrnb);
1669         }
1670     }
1671
1672     if (stepWork.computeNonbondedForces && !useOrEmulateGpuNb)
1673     {
1674         if (havePPDomainDecomposition(cr))
1675         {
1676             do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1677                          nrnb, wcycle);
1678         }
1679
1680         if (stepWork.computeForces)
1681         {
1682             /* Add all the non-bonded force to the normal force array.
1683              * This can be split into a local and a non-local part when overlapping
1684              * communication with calculation with domain decomposition.
1685              */
1686             wallcycle_stop(wcycle, ewcFORCE);
1687             nbv->atomdata_add_nbat_f_to_f(AtomLocality::All,
1688                                           forceOutNonbonded->forceWithShiftForces().force());
1689             wallcycle_start_nocount(wcycle, ewcFORCE);
1690         }
1691
1692         /* If there are multiple fshift output buffers we need to reduce them */
1693         if (stepWork.computeVirial)
1694         {
1695             /* This is not in a subcounter because it takes a
1696                negligible and constant-sized amount of time */
1697             nbnxn_atomdata_add_nbat_fshift_to_fshift(
1698                     *nbv->nbat, forceOutNonbonded->forceWithShiftForces().shiftForces());
1699         }
1700     }
1701
1702     // TODO Force flags should include haveFreeEnergyWork for this domain
1703     if (stepWork.useGpuXHalo && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1704     {
1705         /* Wait for non-local coordinate data to be copied from device */
1706         stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1707     }
1708
1709     // Compute wall interactions, when present.
1710     // Note: should be moved to special forces.
1711     if (inputrec->nwall && stepWork.computeNonbondedForces)
1712     {
1713         /* foreign lambda component for walls */
1714         real dvdl_walls = do_walls(*inputrec, *fr, box, *mdatoms, x.unpaddedConstArrayRef(),
1715                                    &forceOutMtsLevel0.forceWithVirial(), lambda[efptVDW],
1716                                    enerd->grpp.ener[egLJSR].data(), nrnb);
1717         enerd->dvdl_lin[efptVDW] += dvdl_walls;
1718     }
1719
1720     if (stepWork.computeListedForces)
1721     {
1722         /* Check whether we need to take into account PBC in listed interactions */
1723         bool needMolPbc = false;
1724         for (const auto& listedForces : fr->listedForces)
1725         {
1726             if (listedForces.haveCpuListedForces(*fr->fcdata))
1727             {
1728                 needMolPbc = fr->bMolPBC;
1729             }
1730         }
1731
1732         t_pbc pbc;
1733
1734         if (needMolPbc)
1735         {
1736             /* Since all atoms are in the rectangular or triclinic unit-cell,
1737              * only single box vector shifts (2 in x) are required.
1738              */
1739             set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1740         }
1741
1742         for (int mtsIndex = 0; mtsIndex < (fr->useMts && stepWork.computeSlowForces ? 2 : 1); mtsIndex++)
1743         {
1744             ListedForces& listedForces = fr->listedForces[mtsIndex];
1745             ForceOutputs& forceOut     = (mtsIndex == 0 ? forceOutMtsLevel0 : *forceOutMtsLevel1);
1746             listedForces.calculate(
1747                     wcycle, box, inputrec->fepvals, cr, ms, x, xWholeMolecules, fr->fcdata.get(),
1748                     hist, &forceOut, fr, &pbc, enerd, nrnb, lambda.data(), mdatoms,
1749                     DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr, stepWork);
1750         }
1751     }
1752
1753     if (stepWork.computeSlowForces)
1754     {
1755         calculateLongRangeNonbondeds(fr, inputrec, cr, nrnb, wcycle, mdatoms,
1756                                      x.unpaddedConstArrayRef(), &forceOutMtsLevel1->forceWithVirial(),
1757                                      enerd, box, lambda.data(), as_rvec_array(dipoleData.muStateAB),
1758                                      stepWork, ddBalanceRegionHandler);
1759     }
1760
1761     wallcycle_stop(wcycle, ewcFORCE);
1762
1763     // VdW dispersion correction, only computed on master rank to avoid double counting
1764     if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1765     {
1766         // Calculate long range corrections to pressure and energy
1767         const DispersionCorrection::Correction correction =
1768                 fr->dispersionCorrection->calculate(box, lambda[efptVDW]);
1769
1770         if (stepWork.computeEnergy)
1771         {
1772             enerd->term[F_DISPCORR] = correction.energy;
1773             enerd->term[F_DVDL_VDW] += correction.dvdl;
1774             enerd->dvdl_lin[efptVDW] += correction.dvdl;
1775         }
1776         if (stepWork.computeVirial)
1777         {
1778             correction.correctVirial(vir_force);
1779             enerd->term[F_PDISPCORR] = correction.pressure;
1780         }
1781     }
1782
1783     computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation, imdSession, pull_work, step, t,
1784                          wcycle, fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda, stepWork,
1785                          &forceOutMtsLevel0.forceWithVirial(), enerd, ed, stepWork.doNeighborSearch);
1786
1787     GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
1788                "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
1789     GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFHalo),
1790                "The schedule below does not allow for nonbonded MTS with GPU halo exchange");
1791     // Will store the amount of cycles spent waiting for the GPU that
1792     // will be later used in the DLB accounting.
1793     float cycles_wait_gpu = 0;
1794     if (useOrEmulateGpuNb && stepWork.computeNonbondedForces)
1795     {
1796         auto& forceWithShiftForces = forceOutNonbonded->forceWithShiftForces();
1797
1798         /* wait for non-local forces (or calculate in emulation mode) */
1799         if (havePPDomainDecomposition(cr))
1800         {
1801             if (simulationWork.useGpuNonbonded)
1802             {
1803                 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1804                         nbv->gpu_nbv, stepWork, AtomLocality::NonLocal, enerd->grpp.ener[egLJSR].data(),
1805                         enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(), wcycle);
1806             }
1807             else
1808             {
1809                 wallcycle_start_nocount(wcycle, ewcFORCE);
1810                 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes,
1811                              step, nrnb, wcycle);
1812                 wallcycle_stop(wcycle, ewcFORCE);
1813             }
1814
1815             if (stepWork.useGpuFBufferOps)
1816             {
1817                 // TODO: move this into DomainLifetimeWorkload, including the second part of the
1818                 // condition The bonded and free energy CPU tasks can have non-local force
1819                 // contributions which are a dependency for the GPU force reduction.
1820                 bool haveNonLocalForceContribInCpuBuffer =
1821                         domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1822
1823                 if (haveNonLocalForceContribInCpuBuffer)
1824                 {
1825                     stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
1826                                               AtomLocality::NonLocal);
1827                 }
1828
1829
1830                 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->execute();
1831
1832                 if (!stepWork.useGpuFHalo)
1833                 {
1834                     // copy from GPU input for dd_move_f()
1835                     stateGpu->copyForcesFromGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
1836                                                 AtomLocality::NonLocal);
1837                 }
1838             }
1839             else
1840             {
1841                 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
1842             }
1843
1844             if (fr->nbv->emulateGpu() && stepWork.computeVirial)
1845             {
1846                 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
1847             }
1848         }
1849     }
1850
1851     /* Combining the forces for multiple time stepping before the halo exchange, when possible,
1852      * avoids an extra halo exchange (when DD is used) and post-processing step.
1853      */
1854     const bool combineMtsForcesBeforeHaloExchange =
1855             (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
1856              && (legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0
1857              && !(stepWork.computeVirial || simulationWork.useGpuNonbonded || useGpuPmeOnThisRank));
1858     if (combineMtsForcesBeforeHaloExchange)
1859     {
1860         const int numAtoms = havePPDomainDecomposition(cr) ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr;
1861         combineMtsForces(numAtoms, force.unpaddedArrayRef(), forceView->forceMtsCombined(),
1862                          inputrec->mtsLevels[1].stepFactor);
1863     }
1864
1865     if (havePPDomainDecomposition(cr))
1866     {
1867         /* We are done with the CPU compute.
1868          * We will now communicate the non-local forces.
1869          * If we use a GPU this will overlap with GPU work, so in that case
1870          * we do not close the DD force balancing region here.
1871          */
1872         ddBalanceRegionHandler.closeAfterForceComputationCpu();
1873
1874         if (stepWork.computeForces)
1875         {
1876
1877             if (stepWork.useGpuFHalo)
1878             {
1879                 if (domainWork.haveCpuLocalForceWork)
1880                 {
1881                     stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
1882                                               AtomLocality::Local);
1883                 }
1884                 communicateGpuHaloForces(*cr, domainWork.haveCpuLocalForceWork);
1885             }
1886             else
1887             {
1888                 if (stepWork.useGpuFBufferOps)
1889                 {
1890                     stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
1891                 }
1892
1893                 // Without MTS or with MTS at slow steps with uncombined forces we need to
1894                 // communicate the fast forces
1895                 if (!fr->useMts || !combineMtsForcesBeforeHaloExchange)
1896                 {
1897                     dd_move_f(cr->dd, &forceOutMtsLevel0.forceWithShiftForces(), wcycle);
1898                 }
1899                 // With MTS we need to communicate the slow or combined (in forceOutMtsLevel1) forces
1900                 if (fr->useMts && stepWork.computeSlowForces)
1901                 {
1902                     dd_move_f(cr->dd, &forceOutMtsLevel1->forceWithShiftForces(), wcycle);
1903                 }
1904             }
1905         }
1906     }
1907
1908     // With both nonbonded and PME offloaded a GPU on the same rank, we use
1909     // an alternating wait/reduction scheme.
1910     bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
1911                              && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
1912     if (alternateGpuWait)
1913     {
1914         alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, forceOutNonbonded,
1915                                     forceOutMtsLevel1, enerd, lambda[efptCOUL], stepWork, wcycle);
1916     }
1917
1918     if (!alternateGpuWait && useGpuPmeOnThisRank)
1919     {
1920         pme_gpu_wait_and_reduce(fr->pmedata, stepWork, wcycle,
1921                                 &forceOutMtsLevel1->forceWithVirial(), enerd, lambda[efptCOUL]);
1922     }
1923
1924     /* Wait for local GPU NB outputs on the non-alternating wait path */
1925     if (!alternateGpuWait && stepWork.computeNonbondedForces && simulationWork.useGpuNonbonded)
1926     {
1927         /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1928          * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1929          * but even with a step of 0.1 ms the difference is less than 1%
1930          * of the step time.
1931          */
1932         const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
1933         const float waitCycles               = Nbnxm::gpu_wait_finish_task(
1934                 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
1935                 enerd->grpp.ener[egCOULSR].data(),
1936                 forceOutNonbonded->forceWithShiftForces().shiftForces(), wcycle);
1937
1938         if (ddBalanceRegionHandler.useBalancingRegion())
1939         {
1940             DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1941             if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
1942             {
1943                 /* We measured few cycles, it could be that the kernel
1944                  * and transfer finished earlier and there was no actual
1945                  * wait time, only API call overhead.
1946                  * Then the actual time could be anywhere between 0 and
1947                  * cycles_wait_est. We will use half of cycles_wait_est.
1948                  */
1949                 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1950             }
1951             ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1952         }
1953     }
1954
1955     if (fr->nbv->emulateGpu())
1956     {
1957         // NOTE: emulation kernel is not included in the balancing region,
1958         // but emulation mode does not target performance anyway
1959         wallcycle_start_nocount(wcycle, ewcFORCE);
1960         do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local,
1961                      DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes, step, nrnb, wcycle);
1962         wallcycle_stop(wcycle, ewcFORCE);
1963     }
1964
1965     // If on GPU PME-PP comms or GPU update path, receive forces from PME before GPU buffer ops
1966     // TODO refactor this and unify with below default-path call to the same function
1967     if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces
1968         && (simulationWork.useGpuPmePpCommunication || simulationWork.useGpuUpdate))
1969     {
1970         /* In case of node-splitting, the PP nodes receive the long-range
1971          * forces, virial and energy from the PME nodes here.
1972          */
1973         pme_receive_force_ener(fr, cr, &forceOutMtsLevel1->forceWithVirial(), enerd,
1974                                simulationWork.useGpuPmePpCommunication,
1975                                stepWork.useGpuPmeFReduction, wcycle);
1976     }
1977
1978
1979     /* Do the nonbonded GPU (or emulation) force buffer reduction
1980      * on the non-alternating path. */
1981     GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
1982                "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
1983     if (useOrEmulateGpuNb && !alternateGpuWait)
1984     {
1985         if (stepWork.useGpuFBufferOps)
1986         {
1987             ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
1988
1989             // Flag to specify whether the CPU force buffer has contributions to
1990             // local atoms. This depends on whether there are CPU-based force tasks
1991             // or when DD is active the halo exchange has resulted in contributions
1992             // from the non-local part.
1993             const bool haveLocalForceContribInCpuBuffer =
1994                     (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
1995
1996             // TODO: move these steps as early as possible:
1997             // - CPU f H2D should be as soon as all CPU-side forces are done
1998             // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
1999             //   before the next CPU task that consumes the forces: vsite spread or update)
2000             // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
2001             //   of the halo exchange. In that case the copy is instead performed above, before the exchange.
2002             //   These should be unified.
2003             if (haveLocalForceContribInCpuBuffer && !stepWork.useGpuFHalo)
2004             {
2005                 // Note: AtomLocality::All is used for the non-DD case because, as in this
2006                 // case copyForcesToGpu() uses a separate stream, it allows overlap of
2007                 // CPU force H2D with GPU force tasks on all streams including those in the
2008                 // local stream which would otherwise be implicit dependencies for the
2009                 // transfer and would not overlap.
2010                 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
2011
2012                 stateGpu->copyForcesToGpu(forceWithShift, locality);
2013             }
2014
2015             if (stepWork.computeNonbondedForces)
2016             {
2017                 fr->gpuForceReduction[gmx::AtomLocality::Local]->execute();
2018             }
2019
2020             // Copy forces to host if they are needed for update or if virtual sites are enabled.
2021             // If there are vsites, we need to copy forces every step to spread vsite forces on host.
2022             // TODO: When the output flags will be included in step workload, this copy can be combined with the
2023             //       copy call done in sim_utils(...) for the output.
2024             // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
2025             //       they should not be copied in do_md(...) for the output.
2026             if (!simulationWork.useGpuUpdate || vsite)
2027             {
2028                 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
2029                 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
2030             }
2031         }
2032         else if (stepWork.computeNonbondedForces)
2033         {
2034             ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2035             nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
2036         }
2037     }
2038
2039     launchGpuEndOfStepTasks(nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork,
2040                             useGpuPmeOnThisRank, step, wcycle);
2041
2042     if (DOMAINDECOMP(cr))
2043     {
2044         dd_force_flop_stop(cr->dd, nrnb);
2045     }
2046
2047     const bool haveCombinedMtsForces = (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2048                                         && combineMtsForcesBeforeHaloExchange);
2049     if (stepWork.computeForces)
2050     {
2051         postProcessForceWithShiftForces(nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutMtsLevel0,
2052                                         vir_force, *mdatoms, *fr, vsite, stepWork);
2053
2054         if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2055         {
2056             postProcessForceWithShiftForces(nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1,
2057                                             vir_force, *mdatoms, *fr, vsite, stepWork);
2058         }
2059     }
2060
2061     // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
2062     if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
2063         && !simulationWork.useGpuUpdate && stepWork.computeSlowForces)
2064     {
2065         /* In case of node-splitting, the PP nodes receive the long-range
2066          * forces, virial and energy from the PME nodes here.
2067          */
2068         pme_receive_force_ener(fr, cr, &forceOutMtsLevel1->forceWithVirial(), enerd,
2069                                simulationWork.useGpuPmePpCommunication, false, wcycle);
2070     }
2071
2072     if (stepWork.computeForces)
2073     {
2074         /* If we don't use MTS or if we already combined the MTS forces before, we only
2075          * need to post-process one ForceOutputs object here, called forceOutCombined,
2076          * otherwise we have to post-process two outputs and then combine them.
2077          */
2078         ForceOutputs& forceOutCombined = (haveCombinedMtsForces ? forceOutMts.value() : forceOutMtsLevel0);
2079         postProcessForces(cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutCombined,
2080                           vir_force, mdatoms, fr, vsite, stepWork);
2081
2082         if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2083         {
2084             postProcessForces(cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1,
2085                               vir_force, mdatoms, fr, vsite, stepWork);
2086
2087             combineMtsForces(mdatoms->homenr, force.unpaddedArrayRef(),
2088                              forceView->forceMtsCombined(), inputrec->mtsLevels[1].stepFactor);
2089         }
2090     }
2091
2092     if (stepWork.computeEnergy)
2093     {
2094         /* Compute the final potential energy terms */
2095         accumulatePotentialEnergies(enerd, lambda, inputrec->fepvals);
2096
2097         if (!EI_TPI(inputrec->eI))
2098         {
2099             checkPotentialEnergyValidity(step, *enerd, *inputrec);
2100         }
2101     }
2102
2103     /* In case we don't have constraints and are using GPUs, the next balancing
2104      * region starts here.
2105      * Some "special" work at the end of do_force_cuts?, such as vsite spread,
2106      * virial calculation and COM pulling, is not thus not included in
2107      * the balance timing, which is ok as most tasks do communication.
2108      */
2109     ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);
2110 }