Pass ArrayRefs to dispatchNonbondedKernel
[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.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
970     flags.computeNonbondedForces =
971             ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded
972             && !(simulationWork.computeNonbondedAtMtsLevel1 && !computeSlowForces);
973     flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
974
975     if (simulationWork.useGpuBufferOps)
976     {
977         GMX_ASSERT(simulationWork.useGpuNonbonded,
978                    "Can only offload buffer ops if nonbonded computation is also offloaded");
979     }
980     flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
981     // on virial steps the CPU reduction path is taken
982     flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
983     flags.useGpuPmeFReduction = flags.computeSlowForces && flags.useGpuFBufferOps && simulationWork.useGpuPme
984                                 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication);
985     flags.useGpuXHalo = simulationWork.useGpuHaloExchange;
986     flags.useGpuFHalo = simulationWork.useGpuHaloExchange && flags.useGpuFBufferOps;
987
988     return flags;
989 }
990
991
992 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
993  *
994  * TODO: eliminate \p useGpuPmeOnThisRank when this is
995  * incorporated in DomainLifetimeWorkload.
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                                     bool                              useGpuPmeOnThisRank,
1003                                     int64_t                           step,
1004                                     gmx_wallcycle*                    wcycle)
1005 {
1006     if (runScheduleWork.simulationWork.useGpuNonbonded && runScheduleWork.stepWork.computeNonbondedForces)
1007     {
1008         /* Launch pruning before buffer clearing because the API overhead of the
1009          * clear kernel launches can leave the GPU idle while it could be running
1010          * the prune kernel.
1011          */
1012         if (nbv->isDynamicPruningStepGpu(step))
1013         {
1014             nbv->dispatchPruneKernelGpu(step);
1015         }
1016
1017         /* now clear the GPU outputs while we finish the step on the CPU */
1018         wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1019         wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1020         Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
1021         wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1022         wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1023     }
1024
1025     if (useGpuPmeOnThisRank)
1026     {
1027         pme_gpu_reinit_computation(pmedata, wcycle);
1028     }
1029
1030     if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
1031     {
1032         // in principle this should be included in the DD balancing region,
1033         // but generally it is infrequent so we'll omit it for the sake of
1034         // simpler code
1035         gpuBonded->waitAccumulateEnergyTerms(enerd);
1036
1037         gpuBonded->clearEnergies();
1038     }
1039 }
1040
1041 //! \brief Data structure to hold dipole-related data and staging arrays
1042 struct DipoleData
1043 {
1044     //! Dipole staging for fast summing over MPI
1045     gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
1046     //! Dipole staging for states A and B (index 0 and 1 resp.)
1047     gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
1048 };
1049
1050
1051 static void reduceAndUpdateMuTot(DipoleData*                   dipoleData,
1052                                  const t_commrec*              cr,
1053                                  const bool                    haveFreeEnergy,
1054                                  gmx::ArrayRef<const real>     lambda,
1055                                  rvec                          muTotal,
1056                                  const DDBalanceRegionHandler& ddBalanceRegionHandler)
1057 {
1058     if (PAR(cr))
1059     {
1060         gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
1061         ddBalanceRegionHandler.reopenRegionCpu();
1062     }
1063     for (int i = 0; i < 2; i++)
1064     {
1065         for (int j = 0; j < DIM; j++)
1066         {
1067             dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
1068         }
1069     }
1070
1071     if (!haveFreeEnergy)
1072     {
1073         copy_rvec(dipoleData->muStateAB[0], muTotal);
1074     }
1075     else
1076     {
1077         for (int j = 0; j < DIM; j++)
1078         {
1079             muTotal[j] = (1.0 - lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)])
1080                                  * dipoleData->muStateAB[0][j]
1081                          + lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]
1082                                    * dipoleData->muStateAB[1][j];
1083         }
1084     }
1085 }
1086
1087 /*! \brief Combines MTS level0 and level1 force buffes into a full and MTS-combined force buffer.
1088  *
1089  * \param[in]     numAtoms        The number of atoms to combine forces for
1090  * \param[in,out] forceMtsLevel0  Input: F_level0, output: F_level0 + F_level1
1091  * \param[in,out] forceMts        Input: F_level1, output: F_level0 + mtsFactor * F_level1
1092  * \param[in]     mtsFactor       The factor between the level0 and level1 time step
1093  */
1094 static void combineMtsForces(const int      numAtoms,
1095                              ArrayRef<RVec> forceMtsLevel0,
1096                              ArrayRef<RVec> forceMts,
1097                              const real     mtsFactor)
1098 {
1099     const int gmx_unused numThreads = gmx_omp_nthreads_get(ModuleMultiThread::Default);
1100 #pragma omp parallel for num_threads(numThreads) schedule(static)
1101     for (int i = 0; i < numAtoms; i++)
1102     {
1103         const RVec forceMtsLevel0Tmp = forceMtsLevel0[i];
1104         forceMtsLevel0[i] += forceMts[i];
1105         forceMts[i] = forceMtsLevel0Tmp + mtsFactor * forceMts[i];
1106     }
1107 }
1108
1109 /*! \brief Setup for the local and non-local GPU force reductions:
1110  * reinitialization plus the registration of forces and dependencies.
1111  *
1112  * \param [in] runScheduleWork               Schedule workload flag structure
1113  * \param [in] cr                            Communication record object
1114  * \param [in] fr                            Force record object
1115  */
1116 static void setupGpuForceReductions(gmx::MdrunScheduleWorkload* runScheduleWork,
1117                                     const t_commrec*            cr,
1118                                     t_forcerec*                 fr)
1119 {
1120
1121     nonbonded_verlet_t*          nbv      = fr->nbv.get();
1122     gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1123
1124     // (re-)initialize local GPU force reduction
1125     const bool accumulate =
1126             runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr);
1127     const int atomStart = 0;
1128     fr->gpuForceReduction[gmx::AtomLocality::Local]->reinit(stateGpu->getForces(),
1129                                                             nbv->getNumAtoms(AtomLocality::Local),
1130                                                             nbv->getGridIndices(),
1131                                                             atomStart,
1132                                                             accumulate,
1133                                                             stateGpu->fReducedOnDevice());
1134
1135     // register forces and add dependencies
1136     fr->gpuForceReduction[gmx::AtomLocality::Local]->registerNbnxmForce(nbv->getGpuForces());
1137
1138     if (runScheduleWork->simulationWork.useGpuPme
1139         && (thisRankHasDuty(cr, DUTY_PME) || runScheduleWork->simulationWork.useGpuPmePpCommunication))
1140     {
1141         DeviceBuffer<gmx::RVec> forcePtr =
1142                 thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1143                                               :                    // PME force buffer on same GPU
1144                         fr->pmePpCommGpu->getGpuForceStagingPtr(); // buffer received from other GPU
1145         fr->gpuForceReduction[gmx::AtomLocality::Local]->registerRvecForce(forcePtr);
1146
1147         GpuEventSynchronizer* const pmeSynchronizer =
1148                 (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1149                                                : // PME force buffer on same GPU
1150                          fr->pmePpCommGpu->getForcesReadySynchronizer()); // buffer received from other GPU
1151
1152         if (GMX_THREAD_MPI)
1153         {
1154             GMX_ASSERT(pmeSynchronizer != nullptr, "PME force ready cuda event should not be NULL");
1155             fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(pmeSynchronizer);
1156         }
1157     }
1158
1159     if ((runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr))
1160         && !runScheduleWork->simulationWork.useGpuHaloExchange)
1161     {
1162         auto forcesReadyLocality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
1163         const bool useGpuForceBufferOps = true;
1164         fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1165                 stateGpu->getForcesReadyOnDeviceEvent(forcesReadyLocality, useGpuForceBufferOps));
1166     }
1167
1168     if (runScheduleWork->simulationWork.useGpuHaloExchange)
1169     {
1170         fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1171                 cr->dd->gpuHaloExchange[0][0]->getForcesReadyOnDeviceEvent());
1172     }
1173
1174     if (havePPDomainDecomposition(cr))
1175     {
1176         // (re-)initialize non-local GPU force reduction
1177         const bool accumulate = runScheduleWork->domainWork.haveCpuBondedWork
1178                                 || runScheduleWork->domainWork.haveFreeEnergyWork;
1179         const int atomStart = dd_numHomeAtoms(*cr->dd);
1180         fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->reinit(stateGpu->getForces(),
1181                                                                    nbv->getNumAtoms(AtomLocality::NonLocal),
1182                                                                    nbv->getGridIndices(),
1183                                                                    atomStart,
1184                                                                    accumulate);
1185
1186         // register forces and add dependencies
1187         fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->registerNbnxmForce(nbv->getGpuForces());
1188         if (runScheduleWork->domainWork.haveCpuBondedWork || runScheduleWork->domainWork.haveFreeEnergyWork)
1189         {
1190             fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->addDependency(
1191                     stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::NonLocal, true));
1192         }
1193     }
1194 }
1195
1196
1197 void do_force(FILE*                               fplog,
1198               const t_commrec*                    cr,
1199               const gmx_multisim_t*               ms,
1200               const t_inputrec&                   inputrec,
1201               gmx::Awh*                           awh,
1202               gmx_enfrot*                         enforcedRotation,
1203               gmx::ImdSession*                    imdSession,
1204               pull_t*                             pull_work,
1205               int64_t                             step,
1206               t_nrnb*                             nrnb,
1207               gmx_wallcycle*                      wcycle,
1208               const gmx_localtop_t*               top,
1209               const matrix                        box,
1210               gmx::ArrayRefWithPadding<gmx::RVec> x,
1211               const history_t*                    hist,
1212               gmx::ForceBuffersView*              forceView,
1213               tensor                              vir_force,
1214               const t_mdatoms*                    mdatoms,
1215               gmx_enerdata_t*                     enerd,
1216               gmx::ArrayRef<const real>           lambda,
1217               t_forcerec*                         fr,
1218               gmx::MdrunScheduleWorkload*         runScheduleWork,
1219               gmx::VirtualSitesHandler*           vsite,
1220               rvec                                muTotal,
1221               double                              t,
1222               gmx_edsam*                          ed,
1223               int                                 legacyFlags,
1224               const DDBalanceRegionHandler&       ddBalanceRegionHandler)
1225 {
1226     auto force = forceView->forceWithPadding();
1227     GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1228                "The size of the force buffer should be at least the number of atoms to compute "
1229                "forces for");
1230
1231     nonbonded_verlet_t*  nbv = fr->nbv.get();
1232     interaction_const_t* ic  = fr->ic.get();
1233
1234     gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1235
1236     const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1237
1238     runScheduleWork->stepWork = setupStepWorkload(
1239             legacyFlags, inputrec.mtsLevels, step, simulationWork, thisRankHasDuty(cr, DUTY_PME));
1240     const StepWorkload& stepWork = runScheduleWork->stepWork;
1241
1242     const bool useGpuPmeOnThisRank =
1243             simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces;
1244
1245     /* At a search step we need to start the first balancing region
1246      * somewhere early inside the step after communication during domain
1247      * decomposition (and not during the previous step as usual).
1248      */
1249     if (stepWork.doNeighborSearch)
1250     {
1251         ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1252     }
1253
1254     clear_mat(vir_force);
1255
1256     if (fr->pbcType != PbcType::No)
1257     {
1258         /* Compute shift vectors every step,
1259          * because of pressure coupling or box deformation!
1260          */
1261         if (stepWork.haveDynamicBox && stepWork.stateChanged)
1262         {
1263             calc_shifts(box, fr->shift_vec);
1264         }
1265
1266         const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1267         const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1268         if (calcCGCM)
1269         {
1270             put_atoms_in_box_omp(fr->pbcType,
1271                                  box,
1272                                  x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1273                                  gmx_omp_nthreads_get(ModuleMultiThread::Default));
1274             inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1275         }
1276     }
1277
1278     nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1279
1280     const bool pmeSendCoordinatesFromGpu =
1281             GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1282     const bool reinitGpuPmePpComms =
1283             GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1284
1285     auto* localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1286                                         ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1287                                                   AtomLocality::Local, simulationWork, stepWork)
1288                                         : nullptr;
1289
1290     // Copy coordinate from the GPU if update is on the GPU and there
1291     // are forces to be computed on the CPU, or for the computation of
1292     // virial, or if host-side data will be transferred from this task
1293     // to a remote task for halo exchange or PME-PP communication. At
1294     // search steps the current coordinates are already on the host,
1295     // hence copy is not needed.
1296     const bool haveHostPmePpComms =
1297             !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1298
1299     GMX_ASSERT(simulationWork.useGpuHaloExchange
1300                        == ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange[0].empty())),
1301                "The GPU halo exchange is active, but it has not been constructed.");
1302     const bool haveHostHaloExchangeComms =
1303             havePPDomainDecomposition(cr) && !simulationWork.useGpuHaloExchange;
1304
1305     bool gmx_used_in_debug haveCopiedXFromGpu = false;
1306     if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1307         && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1308             || haveHostPmePpComms || haveHostHaloExchangeComms || simulationWork.computeMuTot))
1309     {
1310         stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1311         haveCopiedXFromGpu = true;
1312     }
1313
1314     // Coordinates on the device are needed if PME or BufferOps are offloaded.
1315     // The local coordinates can be copied right away.
1316     // NOTE: Consider moving this copy to right after they are updated and constrained,
1317     //       if the later is not offloaded.
1318     if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1319     {
1320         if (stepWork.doNeighborSearch)
1321         {
1322             // TODO refactor this to do_md, after partitioning.
1323             stateGpu->reinit(mdatoms->homenr,
1324                              cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1325             if (useGpuPmeOnThisRank)
1326             {
1327                 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1328                 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1329             }
1330         }
1331         // We need to copy coordinates when:
1332         // 1. Update is not offloaded
1333         // 2. The buffers were reinitialized on search step
1334         if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1335         {
1336             GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1337             stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1338         }
1339     }
1340
1341     if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces)
1342     {
1343         /* Send particle coordinates to the pme nodes */
1344         if (!pmeSendCoordinatesFromGpu && !stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1345         {
1346             GMX_ASSERT(haveCopiedXFromGpu,
1347                        "a wait should only be triggered if copy has been scheduled");
1348             stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1349         }
1350
1351         gmx_pme_send_coordinates(fr,
1352                                  cr,
1353                                  box,
1354                                  as_rvec_array(x.unpaddedArrayRef().data()),
1355                                  lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1356                                  lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1357                                  (stepWork.computeVirial || stepWork.computeEnergy),
1358                                  step,
1359                                  simulationWork.useGpuPmePpCommunication,
1360                                  reinitGpuPmePpComms,
1361                                  pmeSendCoordinatesFromGpu,
1362                                  localXReadyOnDevice,
1363                                  wcycle);
1364     }
1365
1366     if (useGpuPmeOnThisRank)
1367     {
1368         launchPmeGpuSpread(fr->pmedata,
1369                            box,
1370                            stepWork,
1371                            localXReadyOnDevice,
1372                            lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1373                            wcycle);
1374     }
1375
1376     const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1377
1378     /* do gridding for pair search */
1379     if (stepWork.doNeighborSearch)
1380     {
1381         if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1382         {
1383             fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1384         }
1385
1386         wallcycle_start(wcycle, WallCycleCounter::NS);
1387         if (!DOMAINDECOMP(cr))
1388         {
1389             const rvec vzero       = { 0.0_real, 0.0_real, 0.0_real };
1390             const rvec boxDiagonal = { box[XX][XX], box[YY][YY], box[ZZ][ZZ] };
1391             wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSGridLocal);
1392             nbnxn_put_on_grid(nbv,
1393                               box,
1394                               0,
1395                               vzero,
1396                               boxDiagonal,
1397                               nullptr,
1398                               { 0, mdatoms->homenr },
1399                               -1,
1400                               fr->cginfo,
1401                               x.unpaddedArrayRef(),
1402                               0,
1403                               nullptr);
1404             wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSGridLocal);
1405         }
1406         else
1407         {
1408             wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSGridNonLocal);
1409             nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1410             wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSGridNonLocal);
1411         }
1412
1413         nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1414                                gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1415                                fr->cginfo);
1416
1417         wallcycle_stop(wcycle, WallCycleCounter::NS);
1418
1419         /* initialize the GPU nbnxm atom data and bonded data structures */
1420         if (simulationWork.useGpuNonbonded)
1421         {
1422             // Note: cycle counting only nononbondeds, gpuBonded counts internally
1423             wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1424             wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1425             Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1426             wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1427             wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1428
1429             if (fr->gpuBonded)
1430             {
1431                 /* Now we put all atoms on the grid, we can assign bonded
1432                  * interactions to the GPU, where the grid order is
1433                  * needed. Also the xq, f and fshift device buffers have
1434                  * been reallocated if needed, so the bonded code can
1435                  * learn about them. */
1436                 // TODO the xq, f, and fshift buffers are now shared
1437                 // resources, so they should be maintained by a
1438                 // higher-level object than the nb module.
1439                 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1440                                                                       top->idef,
1441                                                                       Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1442                                                                       Nbnxm::gpu_get_f(nbv->gpu_nbv),
1443                                                                       Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1444             }
1445         }
1446
1447         // Need to run after the GPU-offload bonded interaction lists
1448         // are set up to be able to determine whether there is bonded work.
1449         runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1450                 inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1451
1452         wallcycle_start_nocount(wcycle, WallCycleCounter::NS);
1453         wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSSearchLocal);
1454         /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1455         nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1456
1457         nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1458
1459         wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSSearchLocal);
1460         wallcycle_stop(wcycle, WallCycleCounter::NS);
1461
1462         if (stepWork.useGpuXBufferOps)
1463         {
1464             nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1465         }
1466
1467         if (simulationWork.useGpuBufferOps)
1468         {
1469             setupGpuForceReductions(runScheduleWork, cr, fr);
1470         }
1471     }
1472     else if (!EI_TPI(inputrec.eI) && stepWork.computeNonbondedForces)
1473     {
1474         if (stepWork.useGpuXBufferOps)
1475         {
1476             GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1477             nbv->convertCoordinatesGpu(AtomLocality::Local, stateGpu->getCoordinates(), localXReadyOnDevice);
1478         }
1479         else
1480         {
1481             if (simulationWork.useGpuUpdate)
1482             {
1483                 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1484                 GMX_ASSERT(haveCopiedXFromGpu,
1485                            "a wait should only be triggered if copy has been scheduled");
1486                 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1487             }
1488             nbv->convertCoordinates(AtomLocality::Local, x.unpaddedArrayRef());
1489         }
1490     }
1491
1492     if (simulationWork.useGpuNonbonded && (stepWork.computeNonbondedForces || domainWork.haveGpuBondedWork))
1493     {
1494         ddBalanceRegionHandler.openBeforeForceComputationGpu();
1495
1496         wallcycle_start(wcycle, WallCycleCounter::LaunchGpu);
1497         wallcycle_sub_start(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1498         Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1499         if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1500         {
1501             Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1502         }
1503         wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1504         wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1505         // with X buffer ops offloaded to the GPU on all but the search steps
1506
1507         // bonded work not split into separate local and non-local, so with DD
1508         // we can only launch the kernel after non-local coordinates have been received.
1509         if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1510         {
1511             fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1512         }
1513
1514         /* launch local nonbonded work on GPU */
1515         wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1516         wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1517         do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1518         wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1519         wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1520     }
1521
1522     if (useGpuPmeOnThisRank)
1523     {
1524         // In PME GPU and mixed mode we launch FFT / gather after the
1525         // X copy/transform to allow overlap as well as after the GPU NB
1526         // launch to avoid FFT launch overhead hijacking the CPU and delaying
1527         // the nonbonded kernel.
1528         launchPmeGpuFftAndGather(fr->pmedata,
1529                                  lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1530                                  wcycle,
1531                                  stepWork);
1532     }
1533
1534     /* Communicate coordinates and sum dipole if necessary +
1535        do non-local pair search */
1536     if (havePPDomainDecomposition(cr))
1537     {
1538         if (stepWork.doNeighborSearch)
1539         {
1540             // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1541             wallcycle_start_nocount(wcycle, WallCycleCounter::NS);
1542             wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSSearchNonLocal);
1543             /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1544             nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1545
1546             nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1547             wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSSearchNonLocal);
1548             wallcycle_stop(wcycle, WallCycleCounter::NS);
1549             // TODO refactor this GPU halo exchange re-initialisation
1550             // to location in do_md where GPU halo exchange is
1551             // constructed at partitioning, after above stateGpu
1552             // re-initialization has similarly been refactored
1553             if (simulationWork.useGpuHaloExchange)
1554             {
1555                 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1556             }
1557         }
1558         else
1559         {
1560             if (stepWork.useGpuXHalo)
1561             {
1562                 // The following must be called after local setCoordinates (which records an event
1563                 // when the coordinate data has been copied to the device).
1564                 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1565
1566                 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1567                 {
1568                     // non-local part of coordinate buffer must be copied back to host for CPU work
1569                     stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1570                 }
1571             }
1572             else
1573             {
1574                 if (simulationWork.useGpuUpdate)
1575                 {
1576                     GMX_ASSERT(haveCopiedXFromGpu,
1577                                "a wait should only be triggered if copy has been scheduled");
1578                     stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1579                 }
1580                 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1581             }
1582
1583             if (stepWork.useGpuXBufferOps)
1584             {
1585                 if (!useGpuPmeOnThisRank && !stepWork.useGpuXHalo)
1586                 {
1587                     stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1588                 }
1589                 nbv->convertCoordinatesGpu(AtomLocality::NonLocal,
1590                                            stateGpu->getCoordinates(),
1591                                            stateGpu->getCoordinatesReadyOnDeviceEvent(
1592                                                    AtomLocality::NonLocal, simulationWork, stepWork));
1593             }
1594             else
1595             {
1596                 nbv->convertCoordinates(AtomLocality::NonLocal, x.unpaddedArrayRef());
1597             }
1598         }
1599
1600         if (simulationWork.useGpuNonbonded)
1601         {
1602
1603             if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1604             {
1605                 wallcycle_start(wcycle, WallCycleCounter::LaunchGpu);
1606                 wallcycle_sub_start(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1607                 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1608                 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1609                 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1610             }
1611
1612             if (domainWork.haveGpuBondedWork)
1613             {
1614                 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1615             }
1616
1617             /* launch non-local nonbonded tasks on GPU */
1618             wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1619             wallcycle_sub_start(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1620             do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1621             wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1622             wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1623         }
1624     }
1625
1626     if (simulationWork.useGpuNonbonded && stepWork.computeNonbondedForces)
1627     {
1628         /* launch D2H copy-back F */
1629         wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1630         wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1631
1632         if (havePPDomainDecomposition(cr))
1633         {
1634             Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1635         }
1636         Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1637         wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1638
1639         if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1640         {
1641             fr->gpuBonded->launchEnergyTransfer();
1642         }
1643         wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1644     }
1645
1646     gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1647     if (fr->wholeMoleculeTransform)
1648     {
1649         xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1650     }
1651
1652     DipoleData dipoleData;
1653
1654     if (simulationWork.computeMuTot)
1655     {
1656         const int start = 0;
1657
1658         if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch)
1659         {
1660             GMX_ASSERT(haveCopiedXFromGpu,
1661                        "a wait should only be triggered if copy has been scheduled");
1662             stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1663         }
1664
1665         /* Calculate total (local) dipole moment in a temporary common array.
1666          * This makes it possible to sum them over nodes faster.
1667          */
1668         gmx::ArrayRef<const gmx::RVec> xRef =
1669                 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1670         calc_mu(start,
1671                 mdatoms->homenr,
1672                 xRef,
1673                 gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1674                 gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
1675                 mdatoms->nChargePerturbed != 0,
1676                 dipoleData.muStaging[0],
1677                 dipoleData.muStaging[1]);
1678
1679         reduceAndUpdateMuTot(
1680                 &dipoleData, cr, (fr->efep != FreeEnergyPerturbationType::No), lambda, muTotal, ddBalanceRegionHandler);
1681     }
1682
1683     /* Reset energies */
1684     reset_enerdata(enerd);
1685
1686     if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1687     {
1688         wallcycle_start(wcycle, WallCycleCounter::PpDuringPme);
1689         dd_force_flop_start(cr->dd, nrnb);
1690     }
1691
1692     // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1693     // this wait ensures that the D2H transfer is complete.
1694     if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1695         && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1696     {
1697         GMX_ASSERT(haveCopiedXFromGpu, "a wait should only be triggered if copy has been scheduled");
1698         stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1699     }
1700
1701     if (inputrec.bRot)
1702     {
1703         wallcycle_start(wcycle, WallCycleCounter::Rot);
1704         do_rotation(cr, enforcedRotation, box, x.unpaddedConstArrayRef(), t, step, stepWork.doNeighborSearch);
1705         wallcycle_stop(wcycle, WallCycleCounter::Rot);
1706     }
1707
1708     /* Start the force cycle counter.
1709      * Note that a different counter is used for dynamic load balancing.
1710      */
1711     wallcycle_start(wcycle, WallCycleCounter::Force);
1712
1713     /* Set up and clear force outputs:
1714      * forceOutMtsLevel0:  everything except what is in the other two outputs
1715      * forceOutMtsLevel1:  PME-mesh and listed-forces group 1
1716      * forceOutNonbonded: non-bonded forces
1717      * Without multiple time stepping all point to the same object.
1718      * With multiple time-stepping the use is different for MTS fast (level0 only) and slow steps.
1719      */
1720     ForceOutputs forceOutMtsLevel0 = setupForceOutputs(
1721             &fr->forceHelperBuffers[0], force, domainWork, stepWork, havePPDomainDecomposition(cr), wcycle);
1722
1723     // Force output for MTS combined forces, only set at level1 MTS steps
1724     std::optional<ForceOutputs> forceOutMts =
1725             (fr->useMts && stepWork.computeSlowForces)
1726                     ? std::optional(setupForceOutputs(&fr->forceHelperBuffers[1],
1727                                                       forceView->forceMtsCombinedWithPadding(),
1728                                                       domainWork,
1729                                                       stepWork,
1730                                                       havePPDomainDecomposition(cr),
1731                                                       wcycle))
1732                     : std::nullopt;
1733
1734     ForceOutputs* forceOutMtsLevel1 =
1735             fr->useMts ? (stepWork.computeSlowForces ? &forceOutMts.value() : nullptr) : &forceOutMtsLevel0;
1736
1737     const bool nonbondedAtMtsLevel1 = runScheduleWork->simulationWork.computeNonbondedAtMtsLevel1;
1738
1739     ForceOutputs* forceOutNonbonded = nonbondedAtMtsLevel1 ? forceOutMtsLevel1 : &forceOutMtsLevel0;
1740
1741     if (inputrec.bPull && pull_have_constraint(*pull_work))
1742     {
1743         clear_pull_forces(pull_work);
1744     }
1745
1746     /* We calculate the non-bonded forces, when done on the CPU, here.
1747      * We do this before calling do_force_lowlevel, because in that
1748      * function, the listed forces are calculated before PME, which
1749      * does communication.  With this order, non-bonded and listed
1750      * force calculation imbalance can be balanced out by the domain
1751      * decomposition load balancing.
1752      */
1753
1754     const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1755
1756     if (!useOrEmulateGpuNb)
1757     {
1758         do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1759     }
1760
1761     if (fr->efep != FreeEnergyPerturbationType::No && stepWork.computeNonbondedForces)
1762     {
1763         /* Calculate the local and non-local free energy interactions here.
1764          * Happens here on the CPU both with and without GPU.
1765          */
1766         nbv->dispatchFreeEnergyKernel(InteractionLocality::Local,
1767                                       *fr,
1768                                       x.unpaddedArrayRef(),
1769                                       &forceOutNonbonded->forceWithShiftForces(),
1770                                       gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1771                                       gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
1772                                       gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1773                                       gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr),
1774                                       inputrec.fepvals.get(),
1775                                       lambda,
1776                                       enerd,
1777                                       stepWork,
1778                                       nrnb);
1779
1780         if (havePPDomainDecomposition(cr))
1781         {
1782             nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal,
1783                                           *fr,
1784                                           x.unpaddedArrayRef(),
1785                                           &forceOutNonbonded->forceWithShiftForces(),
1786                                           gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1787                                           gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
1788                                           gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1789                                           gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr),
1790                                           inputrec.fepvals.get(),
1791                                           lambda,
1792                                           enerd,
1793                                           stepWork,
1794                                           nrnb);
1795         }
1796     }
1797
1798     if (stepWork.computeNonbondedForces && !useOrEmulateGpuNb)
1799     {
1800         if (havePPDomainDecomposition(cr))
1801         {
1802             do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1803         }
1804
1805         if (stepWork.computeForces)
1806         {
1807             /* Add all the non-bonded force to the normal force array.
1808              * This can be split into a local and a non-local part when overlapping
1809              * communication with calculation with domain decomposition.
1810              */
1811             wallcycle_stop(wcycle, WallCycleCounter::Force);
1812             nbv->atomdata_add_nbat_f_to_f(AtomLocality::All,
1813                                           forceOutNonbonded->forceWithShiftForces().force());
1814             wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
1815         }
1816
1817         /* If there are multiple fshift output buffers we need to reduce them */
1818         if (stepWork.computeVirial)
1819         {
1820             /* This is not in a subcounter because it takes a
1821                negligible and constant-sized amount of time */
1822             nbnxn_atomdata_add_nbat_fshift_to_fshift(
1823                     *nbv->nbat, forceOutNonbonded->forceWithShiftForces().shiftForces());
1824         }
1825     }
1826
1827     // TODO Force flags should include haveFreeEnergyWork for this domain
1828     if (stepWork.useGpuXHalo && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1829     {
1830         wallcycle_stop(wcycle, WallCycleCounter::Force);
1831         /* Wait for non-local coordinate data to be copied from device */
1832         stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1833         wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
1834     }
1835
1836     // Compute wall interactions, when present.
1837     // Note: should be moved to special forces.
1838     if (inputrec.nwall && stepWork.computeNonbondedForces)
1839     {
1840         /* foreign lambda component for walls */
1841         real dvdl_walls = do_walls(inputrec,
1842                                    *fr,
1843                                    box,
1844                                    gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1845                                    gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr),
1846                                    gmx::arrayRefFromArray(mdatoms->cENER, mdatoms->nr),
1847                                    mdatoms->homenr,
1848                                    mdatoms->nPerturbed,
1849                                    x.unpaddedConstArrayRef(),
1850                                    &forceOutMtsLevel0.forceWithVirial(),
1851                                    lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1852                                    enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR],
1853                                    nrnb);
1854         enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += dvdl_walls;
1855     }
1856
1857     if (stepWork.computeListedForces)
1858     {
1859         /* Check whether we need to take into account PBC in listed interactions */
1860         bool needMolPbc = false;
1861         for (const auto& listedForces : fr->listedForces)
1862         {
1863             if (listedForces.haveCpuListedForces(*fr->fcdata))
1864             {
1865                 needMolPbc = fr->bMolPBC;
1866             }
1867         }
1868
1869         t_pbc pbc;
1870
1871         if (needMolPbc)
1872         {
1873             /* Since all atoms are in the rectangular or triclinic unit-cell,
1874              * only single box vector shifts (2 in x) are required.
1875              */
1876             set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1877         }
1878
1879         for (int mtsIndex = 0; mtsIndex < (fr->useMts && stepWork.computeSlowForces ? 2 : 1); mtsIndex++)
1880         {
1881             ListedForces& listedForces = fr->listedForces[mtsIndex];
1882             ForceOutputs& forceOut     = (mtsIndex == 0 ? forceOutMtsLevel0 : *forceOutMtsLevel1);
1883             listedForces.calculate(wcycle,
1884                                    box,
1885                                    inputrec.fepvals.get(),
1886                                    cr,
1887                                    ms,
1888                                    x,
1889                                    xWholeMolecules,
1890                                    fr->fcdata.get(),
1891                                    hist,
1892                                    &forceOut,
1893                                    fr,
1894                                    &pbc,
1895                                    enerd,
1896                                    nrnb,
1897                                    lambda,
1898                                    mdatoms,
1899                                    DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr,
1900                                    stepWork);
1901         }
1902     }
1903
1904     if (stepWork.computeSlowForces)
1905     {
1906         calculateLongRangeNonbondeds(fr,
1907                                      inputrec,
1908                                      cr,
1909                                      nrnb,
1910                                      wcycle,
1911                                      mdatoms,
1912                                      x.unpaddedConstArrayRef(),
1913                                      &forceOutMtsLevel1->forceWithVirial(),
1914                                      enerd,
1915                                      box,
1916                                      lambda,
1917                                      dipoleData.muStateAB,
1918                                      stepWork,
1919                                      ddBalanceRegionHandler);
1920     }
1921
1922     wallcycle_stop(wcycle, WallCycleCounter::Force);
1923
1924     // VdW dispersion correction, only computed on master rank to avoid double counting
1925     if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1926     {
1927         // Calculate long range corrections to pressure and energy
1928         const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
1929                 box, lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)]);
1930
1931         if (stepWork.computeEnergy)
1932         {
1933             enerd->term[F_DISPCORR] = correction.energy;
1934             enerd->term[F_DVDL_VDW] += correction.dvdl;
1935             enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += correction.dvdl;
1936         }
1937         if (stepWork.computeVirial)
1938         {
1939             correction.correctVirial(vir_force);
1940             enerd->term[F_PDISPCORR] = correction.pressure;
1941         }
1942     }
1943
1944     computeSpecialForces(fplog,
1945                          cr,
1946                          inputrec,
1947                          awh,
1948                          enforcedRotation,
1949                          imdSession,
1950                          pull_work,
1951                          step,
1952                          t,
1953                          wcycle,
1954                          fr->forceProviders,
1955                          box,
1956                          x.unpaddedArrayRef(),
1957                          mdatoms,
1958                          lambda,
1959                          stepWork,
1960                          &forceOutMtsLevel0.forceWithVirial(),
1961                          forceOutMtsLevel1 ? &forceOutMtsLevel1->forceWithVirial() : nullptr,
1962                          enerd,
1963                          ed,
1964                          stepWork.doNeighborSearch);
1965
1966     if (havePPDomainDecomposition(cr) && stepWork.computeForces && stepWork.useGpuFHalo
1967         && domainWork.haveCpuLocalForceWork)
1968     {
1969         stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(), AtomLocality::Local);
1970     }
1971
1972     GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
1973                "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
1974     GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFHalo),
1975                "The schedule below does not allow for nonbonded MTS with GPU halo exchange");
1976     // Will store the amount of cycles spent waiting for the GPU that
1977     // will be later used in the DLB accounting.
1978     float cycles_wait_gpu = 0;
1979     if (useOrEmulateGpuNb && stepWork.computeNonbondedForces)
1980     {
1981         auto& forceWithShiftForces = forceOutNonbonded->forceWithShiftForces();
1982
1983         /* wait for non-local forces (or calculate in emulation mode) */
1984         if (havePPDomainDecomposition(cr))
1985         {
1986             if (simulationWork.useGpuNonbonded)
1987             {
1988                 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1989                         nbv->gpu_nbv,
1990                         stepWork,
1991                         AtomLocality::NonLocal,
1992                         enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
1993                         enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
1994                         forceWithShiftForces.shiftForces(),
1995                         wcycle);
1996             }
1997             else
1998             {
1999                 wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
2000                 do_nb_verlet(
2001                         fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes, step, nrnb, wcycle);
2002                 wallcycle_stop(wcycle, WallCycleCounter::Force);
2003             }
2004
2005             if (stepWork.useGpuFBufferOps)
2006             {
2007                 // TODO: move this into DomainLifetimeWorkload, including the second part of the
2008                 // condition The bonded and free energy CPU tasks can have non-local force
2009                 // contributions which are a dependency for the GPU force reduction.
2010                 bool haveNonLocalForceContribInCpuBuffer =
2011                         domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
2012
2013                 if (haveNonLocalForceContribInCpuBuffer)
2014                 {
2015                     stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2016                                               AtomLocality::NonLocal);
2017                 }
2018
2019
2020                 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->execute();
2021
2022                 if (!stepWork.useGpuFHalo)
2023                 {
2024                     // copy from GPU input for dd_move_f()
2025                     stateGpu->copyForcesFromGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2026                                                 AtomLocality::NonLocal);
2027                 }
2028             }
2029             else
2030             {
2031                 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
2032             }
2033
2034             if (fr->nbv->emulateGpu() && stepWork.computeVirial)
2035             {
2036                 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
2037             }
2038         }
2039     }
2040
2041     /* Combining the forces for multiple time stepping before the halo exchange, when possible,
2042      * avoids an extra halo exchange (when DD is used) and post-processing step.
2043      */
2044     const bool combineMtsForcesBeforeHaloExchange =
2045             (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2046              && (legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0
2047              && !(stepWork.computeVirial || simulationWork.useGpuNonbonded || useGpuPmeOnThisRank));
2048     if (combineMtsForcesBeforeHaloExchange)
2049     {
2050         const int numAtoms = havePPDomainDecomposition(cr) ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr;
2051         combineMtsForces(numAtoms,
2052                          force.unpaddedArrayRef(),
2053                          forceView->forceMtsCombined(),
2054                          inputrec.mtsLevels[1].stepFactor);
2055     }
2056
2057     if (havePPDomainDecomposition(cr))
2058     {
2059         /* We are done with the CPU compute.
2060          * We will now communicate the non-local forces.
2061          * If we use a GPU this will overlap with GPU work, so in that case
2062          * we do not close the DD force balancing region here.
2063          */
2064         ddBalanceRegionHandler.closeAfterForceComputationCpu();
2065
2066         if (stepWork.computeForces)
2067         {
2068
2069             if (stepWork.useGpuFHalo)
2070             {
2071                 // If there exist CPU forces, data from halo exchange should accumulate into these
2072                 bool accumulateForces = domainWork.haveCpuLocalForceWork;
2073                 if (!accumulateForces)
2074                 {
2075                     // Force halo exchange will set a subset of local atoms with remote non-local data
2076                     // First clear local portion of force array, so that untouched atoms are zero
2077                     stateGpu->clearForcesOnGpu(AtomLocality::Local);
2078                 }
2079                 communicateGpuHaloForces(*cr, accumulateForces);
2080             }
2081             else
2082             {
2083                 if (stepWork.useGpuFBufferOps)
2084                 {
2085                     stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
2086                 }
2087
2088                 // Without MTS or with MTS at slow steps with uncombined forces we need to
2089                 // communicate the fast forces
2090                 if (!fr->useMts || !combineMtsForcesBeforeHaloExchange)
2091                 {
2092                     dd_move_f(cr->dd, &forceOutMtsLevel0.forceWithShiftForces(), wcycle);
2093                 }
2094                 // With MTS we need to communicate the slow or combined (in forceOutMtsLevel1) forces
2095                 if (fr->useMts && stepWork.computeSlowForces)
2096                 {
2097                     dd_move_f(cr->dd, &forceOutMtsLevel1->forceWithShiftForces(), wcycle);
2098                 }
2099             }
2100         }
2101     }
2102
2103     // With both nonbonded and PME offloaded a GPU on the same rank, we use
2104     // an alternating wait/reduction scheme.
2105     bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
2106                              && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
2107     if (alternateGpuWait)
2108     {
2109         alternatePmeNbGpuWaitReduce(fr->nbv.get(),
2110                                     fr->pmedata,
2111                                     forceOutNonbonded,
2112                                     forceOutMtsLevel1,
2113                                     enerd,
2114                                     lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
2115                                     stepWork,
2116                                     wcycle);
2117     }
2118
2119     if (!alternateGpuWait && useGpuPmeOnThisRank)
2120     {
2121         pme_gpu_wait_and_reduce(fr->pmedata,
2122                                 stepWork,
2123                                 wcycle,
2124                                 &forceOutMtsLevel1->forceWithVirial(),
2125                                 enerd,
2126                                 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]);
2127     }
2128
2129     /* Wait for local GPU NB outputs on the non-alternating wait path */
2130     if (!alternateGpuWait && stepWork.computeNonbondedForces && simulationWork.useGpuNonbonded)
2131     {
2132         /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
2133          * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
2134          * but even with a step of 0.1 ms the difference is less than 1%
2135          * of the step time.
2136          */
2137         const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
2138         const float waitCycles               = Nbnxm::gpu_wait_finish_task(
2139                 nbv->gpu_nbv,
2140                 stepWork,
2141                 AtomLocality::Local,
2142                 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
2143                 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
2144                 forceOutNonbonded->forceWithShiftForces().shiftForces(),
2145                 wcycle);
2146
2147         if (ddBalanceRegionHandler.useBalancingRegion())
2148         {
2149             DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
2150             if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
2151             {
2152                 /* We measured few cycles, it could be that the kernel
2153                  * and transfer finished earlier and there was no actual
2154                  * wait time, only API call overhead.
2155                  * Then the actual time could be anywhere between 0 and
2156                  * cycles_wait_est. We will use half of cycles_wait_est.
2157                  */
2158                 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
2159             }
2160             ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
2161         }
2162     }
2163
2164     if (fr->nbv->emulateGpu())
2165     {
2166         // NOTE: emulation kernel is not included in the balancing region,
2167         // but emulation mode does not target performance anyway
2168         wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
2169         do_nb_verlet(fr,
2170                      ic,
2171                      enerd,
2172                      stepWork,
2173                      InteractionLocality::Local,
2174                      DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
2175                      step,
2176                      nrnb,
2177                      wcycle);
2178         wallcycle_stop(wcycle, WallCycleCounter::Force);
2179     }
2180
2181     // If on GPU PME-PP comms path, receive forces from PME before GPU buffer ops
2182     // TODO refactor this and unify with below default-path call to the same function
2183     if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces
2184         && simulationWork.useGpuPmePpCommunication)
2185     {
2186         /* In case of node-splitting, the PP nodes receive the long-range
2187          * forces, virial and energy from the PME nodes here.
2188          */
2189         pme_receive_force_ener(fr,
2190                                cr,
2191                                &forceOutMtsLevel1->forceWithVirial(),
2192                                enerd,
2193                                simulationWork.useGpuPmePpCommunication,
2194                                stepWork.useGpuPmeFReduction,
2195                                wcycle);
2196     }
2197
2198
2199     /* Do the nonbonded GPU (or emulation) force buffer reduction
2200      * on the non-alternating path. */
2201     GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
2202                "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
2203     if (useOrEmulateGpuNb && !alternateGpuWait)
2204     {
2205         if (stepWork.useGpuFBufferOps)
2206         {
2207             ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2208
2209             // Flag to specify whether the CPU force buffer has contributions to
2210             // local atoms. This depends on whether there are CPU-based force tasks
2211             // or when DD is active the halo exchange has resulted in contributions
2212             // from the non-local part.
2213             const bool haveLocalForceContribInCpuBuffer =
2214                     (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
2215
2216             // TODO: move these steps as early as possible:
2217             // - CPU f H2D should be as soon as all CPU-side forces are done
2218             // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
2219             //   before the next CPU task that consumes the forces: vsite spread or update)
2220             // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
2221             //   of the halo exchange. In that case the copy is instead performed above, before the exchange.
2222             //   These should be unified.
2223             if (haveLocalForceContribInCpuBuffer && !stepWork.useGpuFHalo)
2224             {
2225                 // Note: AtomLocality::All is used for the non-DD case because, as in this
2226                 // case copyForcesToGpu() uses a separate stream, it allows overlap of
2227                 // CPU force H2D with GPU force tasks on all streams including those in the
2228                 // local stream which would otherwise be implicit dependencies for the
2229                 // transfer and would not overlap.
2230                 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
2231
2232                 stateGpu->copyForcesToGpu(forceWithShift, locality);
2233             }
2234
2235             if (stepWork.computeNonbondedForces)
2236             {
2237                 fr->gpuForceReduction[gmx::AtomLocality::Local]->execute();
2238             }
2239
2240             // Copy forces to host if they are needed for update or if virtual sites are enabled.
2241             // If there are vsites, we need to copy forces every step to spread vsite forces on host.
2242             // TODO: When the output flags will be included in step workload, this copy can be combined with the
2243             //       copy call done in sim_utils(...) for the output.
2244             // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
2245             //       they should not be copied in do_md(...) for the output.
2246             if (!simulationWork.useGpuUpdate
2247                 || (simulationWork.useGpuUpdate && DOMAINDECOMP(cr) && haveHostPmePpComms) || vsite)
2248             {
2249                 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
2250                 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
2251             }
2252         }
2253         else if (stepWork.computeNonbondedForces)
2254         {
2255             ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2256             nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
2257         }
2258     }
2259
2260     launchGpuEndOfStepTasks(
2261             nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork, useGpuPmeOnThisRank, step, wcycle);
2262
2263     if (DOMAINDECOMP(cr))
2264     {
2265         dd_force_flop_stop(cr->dd, nrnb);
2266     }
2267
2268     const bool haveCombinedMtsForces = (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2269                                         && combineMtsForcesBeforeHaloExchange);
2270     if (stepWork.computeForces)
2271     {
2272         postProcessForceWithShiftForces(
2273                 nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutMtsLevel0, vir_force, *mdatoms, *fr, vsite, stepWork);
2274
2275         if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2276         {
2277             postProcessForceWithShiftForces(
2278                     nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, *mdatoms, *fr, vsite, stepWork);
2279         }
2280     }
2281
2282     // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
2283     if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
2284         && stepWork.computeSlowForces)
2285     {
2286         /* In case of node-splitting, the PP nodes receive the long-range
2287          * forces, virial and energy from the PME nodes here.
2288          */
2289         pme_receive_force_ener(fr,
2290                                cr,
2291                                &forceOutMtsLevel1->forceWithVirial(),
2292                                enerd,
2293                                simulationWork.useGpuPmePpCommunication,
2294                                false,
2295                                wcycle);
2296     }
2297
2298     if (stepWork.computeForces)
2299     {
2300         /* If we don't use MTS or if we already combined the MTS forces before, we only
2301          * need to post-process one ForceOutputs object here, called forceOutCombined,
2302          * otherwise we have to post-process two outputs and then combine them.
2303          */
2304         ForceOutputs& forceOutCombined = (haveCombinedMtsForces ? forceOutMts.value() : forceOutMtsLevel0);
2305         postProcessForces(
2306                 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutCombined, vir_force, mdatoms, fr, vsite, stepWork);
2307
2308         if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2309         {
2310             postProcessForces(
2311                     cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, mdatoms, fr, vsite, stepWork);
2312
2313             combineMtsForces(mdatoms->homenr,
2314                              force.unpaddedArrayRef(),
2315                              forceView->forceMtsCombined(),
2316                              inputrec.mtsLevels[1].stepFactor);
2317         }
2318     }
2319
2320     if (stepWork.computeEnergy)
2321     {
2322         /* Compute the final potential energy terms */
2323         accumulatePotentialEnergies(enerd, lambda, inputrec.fepvals.get());
2324
2325         if (!EI_TPI(inputrec.eI))
2326         {
2327             checkPotentialEnergyValidity(step, *enerd, inputrec);
2328         }
2329     }
2330
2331     /* In case we don't have constraints and are using GPUs, the next balancing
2332      * region starts here.
2333      * Some "special" work at the end of do_force_cuts?, such as vsite spread,
2334      * virial calculation and COM pulling, is not thus not included in
2335      * the balance timing, which is ok as most tasks do communication.
2336      */
2337     ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);
2338 }