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