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