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