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