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