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