Centralize management of forces ready on GPU event
[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]  useGpuForceReduction True if GPU-based force reduction is active this step
617  * \param[in]  wcycle               The wallcycle structure
618  */
619 static inline void launchPmeGpuSpread(gmx_pme_t          *pmedata,
620                                       const matrix        box,
621                                       const StepWorkload &stepWork,
622                                       int                 pmeFlags,
623                                       bool                useGpuForceReduction,
624                                       gmx_wallcycle_t     wcycle)
625 {
626     pme_gpu_prepare_computation(pmedata, stepWork.haveDynamicBox, box, wcycle, pmeFlags, useGpuForceReduction);
627     pme_gpu_launch_spread(pmedata, wcycle);
628 }
629
630 /*! \brief Launch the FFT and gather stages of PME GPU
631  *
632  * This function only implements setting the output forces (no accumulation).
633  *
634  * \param[in]  pmedata        The PME structure
635  * \param[in]  wcycle         The wallcycle structure
636  */
637 static void launchPmeGpuFftAndGather(gmx_pme_t        *pmedata,
638                                      gmx_wallcycle_t   wcycle)
639 {
640     pme_gpu_launch_complex_transforms(pmedata, wcycle);
641     pme_gpu_launch_gather(pmedata, wcycle, PmeForceOutputHandling::Set);
642 }
643
644 /*! \brief
645  *  Polling wait for either of the PME or nonbonded GPU tasks.
646  *
647  * Instead of a static order in waiting for GPU tasks, this function
648  * polls checking which of the two tasks completes first, and does the
649  * associated force buffer reduction overlapped with the other task.
650  * By doing that, unlike static scheduling order, it can always overlap
651  * one of the reductions, regardless of the GPU task completion order.
652  *
653  * \param[in]     nbv              Nonbonded verlet structure
654  * \param[in,out] pmedata          PME module data
655  * \param[in,out] forceOutputs     Output buffer for the forces and virial
656  * \param[in,out] enerd            Energy data structure results are reduced into
657  * \param[in]     stepWork         Step schedule flags
658  * \param[in]     pmeFlags         PME flags
659  * \param[in]     wcycle           The wallcycle structure
660  */
661 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t *nbv,
662                                         gmx_pme_t          *pmedata,
663                                         gmx::ForceOutputs  *forceOutputs,
664                                         gmx_enerdata_t     *enerd,
665                                         const StepWorkload &stepWork,
666                                         int                 pmeFlags,
667                                         gmx_wallcycle_t     wcycle)
668 {
669     bool isPmeGpuDone = false;
670     bool isNbGpuDone  = false;
671
672
673
674     gmx::ForceWithShiftForces      &forceWithShiftForces = forceOutputs->forceWithShiftForces();
675     gmx::ForceWithVirial           &forceWithVirial      = forceOutputs->forceWithVirial();
676
677     gmx::ArrayRef<const gmx::RVec>  pmeGpuForces;
678
679     while (!isPmeGpuDone || !isNbGpuDone)
680     {
681         if (!isPmeGpuDone)
682         {
683             GpuTaskCompletion completionType = (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
684             isPmeGpuDone = pme_gpu_try_finish_task(pmedata, pmeFlags, wcycle, &forceWithVirial, enerd, completionType);
685         }
686
687         if (!isNbGpuDone)
688         {
689             GpuTaskCompletion completionType = (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
690             isNbGpuDone = Nbnxm::gpu_try_finish_task(nbv->gpu_nbv,
691                                                      stepWork,
692                                                      Nbnxm::AtomLocality::Local,
693                                                      enerd->grpp.ener[egLJSR].data(),
694                                                      enerd->grpp.ener[egCOULSR].data(),
695                                                      forceWithShiftForces.shiftForces(), completionType, wcycle);
696
697             if (isNbGpuDone)
698             {
699                 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
700                                               forceWithShiftForces.force());
701             }
702         }
703     }
704 }
705
706 /*! \brief Set up the different force buffers; also does clearing.
707  *
708  * \param[in] fr        force record pointer
709  * \param[in] pull_work The pull work object.
710  * \param[in] inputrec  input record
711  * \param[in] force     force array
712  * \param[in] stepWork  Step schedule flags
713  * \param[out] wcycle   wallcycle recording structure
714  *
715  * \returns             Cleared force output structure
716  */
717 static ForceOutputs
718 setupForceOutputs(t_forcerec                          *fr,
719                   pull_t                              *pull_work,
720                   const t_inputrec                    &inputrec,
721                   gmx::ArrayRefWithPadding<gmx::RVec>  force,
722                   const StepWorkload                  &stepWork,
723                   gmx_wallcycle_t                      wcycle)
724 {
725     wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
726
727     /* NOTE: We assume fr->shiftForces is all zeros here */
728     gmx::ForceWithShiftForces forceWithShiftForces(force, stepWork.computeVirial, fr->shiftForces);
729
730     if (stepWork.computeForces)
731     {
732         /* Clear the short- and long-range forces */
733         clear_rvecs_omp(fr->natoms_force_constr,
734                         as_rvec_array(forceWithShiftForces.force().data()));
735     }
736
737     /* If we need to compute the virial, we might need a separate
738      * force buffer for algorithms for which the virial is calculated
739      * directly, such as PME. Otherwise, forceWithVirial uses the
740      * the same force (f in legacy calls) buffer as other algorithms.
741      */
742     const bool useSeparateForceWithVirialBuffer = (stepWork.computeForces &&
743                                                    (stepWork.computeVirial && fr->haveDirectVirialContributions));
744     /* forceWithVirial uses the local atom range only */
745     gmx::ForceWithVirial forceWithVirial(useSeparateForceWithVirialBuffer ?
746                                          fr->forceBufferForDirectVirialContributions : force.unpaddedArrayRef(),
747                                          stepWork.computeVirial);
748
749     if (useSeparateForceWithVirialBuffer)
750     {
751         /* TODO: update comment
752          * We only compute forces on local atoms. Note that vsites can
753          * spread to non-local atoms, but that part of the buffer is
754          * cleared separately in the vsite spreading code.
755          */
756         clear_rvecs_omp(forceWithVirial.force_.size(), as_rvec_array(forceWithVirial.force_.data()));
757     }
758
759     if (inputrec.bPull && pull_have_constraint(pull_work))
760     {
761         clear_pull_forces(pull_work);
762     }
763
764     wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
765
766     return ForceOutputs(forceWithShiftForces, forceWithVirial);
767 }
768
769
770 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
771  */
772 static DomainLifetimeWorkload
773 setupDomainLifetimeWorkload(const t_inputrec       &inputrec,
774                             const t_forcerec       &fr,
775                             const pull_t           *pull_work,
776                             const gmx_edsam        *ed,
777                             const t_idef           &idef,
778                             const t_fcdata         &fcd,
779                             const t_mdatoms        &mdatoms,
780                             const StepWorkload     &stepWork)
781 {
782     DomainLifetimeWorkload domainWork;
783     // Note that haveSpecialForces is constant over the whole run
784     domainWork.haveSpecialForces      = haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
785     domainWork.haveCpuBondedWork      = haveCpuBondeds(fr);
786     domainWork.haveGpuBondedWork      = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
787     domainWork.haveRestraintsWork     = havePositionRestraints(idef, fcd);
788     domainWork.haveCpuListedForceWork = haveCpuListedForces(fr, idef, fcd);
789     // Note that haveFreeEnergyWork is constant over the whole run
790     domainWork.haveFreeEnergyWork     = (fr.efep != efepNO && mdatoms.nPerturbed != 0);
791     return domainWork;
792 }
793
794 /*! \brief Set up force flag stuct from the force bitmask.
795  *
796  * \param[in]      legacyFlags          Force bitmask flags used to construct the new flags
797  * \param[in]      isNonbondedOn        Global override, if false forces to turn off all nonbonded calculation.
798  * \returns New Stepworkload description.
799  */
800 static StepWorkload
801 setupStepWorkload(const int     legacyFlags,
802                   const bool    isNonbondedOn)
803 {
804     StepWorkload flags;
805     flags.stateChanged           = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
806     flags.haveDynamicBox         = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
807     flags.doNeighborSearch       = ((legacyFlags & GMX_FORCE_NS) != 0);
808     flags.computeVirial          = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
809     flags.computeEnergy          = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
810     flags.computeForces          = ((legacyFlags & GMX_FORCE_FORCES) != 0);
811     flags.computeListedForces    = ((legacyFlags & GMX_FORCE_LISTED) != 0);
812     flags.computeNonbondedForces = ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && isNonbondedOn;
813     flags.computeDhdl            = ((legacyFlags & GMX_FORCE_DHDL) != 0);
814     return flags;
815 }
816
817
818 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
819  *
820  * TODO: eliminate the \p useGpuNonbonded and \p useGpuNonbonded when these are
821  * incorporated in DomainLifetimeWorkload.
822  */
823 static void
824 launchGpuEndOfStepTasks(nonbonded_verlet_t               *nbv,
825                         gmx::GpuBonded                   *gpuBonded,
826                         gmx_pme_t                        *pmedata,
827                         gmx_enerdata_t                   *enerd,
828                         const gmx::MdrunScheduleWorkload &runScheduleWork,
829                         bool                              useGpuNonbonded,
830                         bool                              useGpuPme,
831                         int64_t                           step,
832                         gmx_wallcycle_t                   wcycle)
833 {
834     if (useGpuNonbonded)
835     {
836         /* Launch pruning before buffer clearing because the API overhead of the
837          * clear kernel launches can leave the GPU idle while it could be running
838          * the prune kernel.
839          */
840         if (nbv->isDynamicPruningStepGpu(step))
841         {
842             nbv->dispatchPruneKernelGpu(step);
843         }
844
845         /* now clear the GPU outputs while we finish the step on the CPU */
846         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
847         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
848         Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
849         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
850         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
851     }
852
853     if (useGpuPme)
854     {
855         pme_gpu_reinit_computation(pmedata, wcycle);
856     }
857
858     if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
859     {
860         // in principle this should be included in the DD balancing region,
861         // but generally it is infrequent so we'll omit it for the sake of
862         // simpler code
863         gpuBonded->waitAccumulateEnergyTerms(enerd);
864
865         gpuBonded->clearEnergies();
866     }
867 }
868
869
870 void do_force(FILE                                     *fplog,
871               const t_commrec                          *cr,
872               const gmx_multisim_t                     *ms,
873               const t_inputrec                         *inputrec,
874               gmx::Awh                                 *awh,
875               gmx_enfrot                               *enforcedRotation,
876               gmx::ImdSession                          *imdSession,
877               pull_t                                   *pull_work,
878               int64_t                                   step,
879               t_nrnb                                   *nrnb,
880               gmx_wallcycle_t                           wcycle,
881               const gmx_localtop_t                     *top,
882               const matrix                              box,
883               gmx::ArrayRefWithPadding<gmx::RVec>       x,
884               history_t                                *hist,
885               gmx::ArrayRefWithPadding<gmx::RVec>       force,
886               tensor                                    vir_force,
887               const t_mdatoms                          *mdatoms,
888               gmx_enerdata_t                           *enerd,
889               t_fcdata                                 *fcd,
890               gmx::ArrayRef<real>                       lambda,
891               t_graph                                  *graph,
892               t_forcerec                               *fr,
893               gmx::MdrunScheduleWorkload               *runScheduleWork,
894               const gmx_vsite_t                        *vsite,
895               rvec                                      mu_tot,
896               double                                    t,
897               gmx_edsam                                *ed,
898               int                                       legacyFlags,
899               const DDBalanceRegionHandler             &ddBalanceRegionHandler)
900 {
901     int                          i, j;
902     double                       mu[2*DIM];
903     nonbonded_verlet_t          *nbv      = fr->nbv.get();
904     interaction_const_t         *ic       = fr->ic;
905     gmx::StatePropagatorDataGpu *stateGpu = fr->stateGpu;
906
907     // TODO remove the code below when the legacy flags are not in use anymore
908     /* modify force flag if not doing nonbonded */
909     if (!fr->bNonbonded)
910     {
911         legacyFlags &= ~GMX_FORCE_NONBONDED;
912     }
913
914     runScheduleWork->stepWork = setupStepWorkload(legacyFlags, fr->bNonbonded);
915     const StepWorkload       &stepWork = runScheduleWork->stepWork;
916
917     const SimulationWorkload &simulationWork = runScheduleWork->simulationWork;
918
919     const bool                useGpuPmeOnThisRank = simulationWork.usePmeGpu && thisRankHasDuty(cr, DUTY_PME);
920     const int                 pmeFlags            = makePmeFlags(stepWork);
921
922     // Switches on whether to use GPU for position and force buffer operations
923     // TODO consider all possible combinations of triggers, and how to combine optimally in each case.
924     const BufferOpsUseGpu useGpuXBufOps = (simulationWork.useGpuBufferOps &&
925                                            simulationWork.useGpuNonbonded && (GMX_GPU == GMX_GPU_CUDA)) ? BufferOpsUseGpu::True : BufferOpsUseGpu::False;;
926     // GPU Force buffer ops are disabled on virial steps, because the virial calc is not yet ported to GPU
927     const BufferOpsUseGpu useGpuFBufOps = ((simulationWork.useGpuBufferOps &&
928                                             simulationWork.useGpuNonbonded && (GMX_GPU == GMX_GPU_CUDA)) &&
929                                            !(stepWork.computeVirial || stepWork.computeEnergy)) ? BufferOpsUseGpu::True : BufferOpsUseGpu::False;
930     // TODO: move / add this flag to the internal PME GPU data structures
931     const bool useGpuPmeFReduction = (useGpuFBufOps == BufferOpsUseGpu::True) &&
932         useGpuPmeOnThisRank; // only supported if this rank is perfoming PME on the GPU
933
934     /* At a search step we need to start the first balancing region
935      * somewhere early inside the step after communication during domain
936      * decomposition (and not during the previous step as usual).
937      */
938     if (stepWork.doNeighborSearch)
939     {
940         ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
941     }
942
943     const int start  = 0;
944     const int homenr = mdatoms->homenr;
945
946     clear_mat(vir_force);
947
948     if (stepWork.stateChanged)
949     {
950         if (inputrecNeedMutot(inputrec))
951         {
952             /* Calculate total (local) dipole moment in a temporary common array.
953              * This makes it possible to sum them over nodes faster.
954              */
955             calc_mu(start, homenr,
956                     x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
957                     mu, mu+DIM);
958         }
959     }
960
961     if (fr->ePBC != epbcNONE)
962     {
963         /* Compute shift vectors every step,
964          * because of pressure coupling or box deformation!
965          */
966         if (stepWork.haveDynamicBox && stepWork.stateChanged)
967         {
968             calc_shifts(box, fr->shift_vec);
969         }
970
971         const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
972         const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
973         if (calcCGCM)
974         {
975             put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr), gmx_omp_nthreads_get(emntDefault));
976             inc_nrnb(nrnb, eNR_SHIFTX, homenr);
977         }
978         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
979         {
980             unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
981         }
982     }
983
984     nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox,
985                                  fr->shift_vec, nbv->nbat.get());
986
987 #if GMX_MPI
988     if (!thisRankHasDuty(cr, DUTY_PME))
989     {
990         /* Send particle coordinates to the pme nodes.
991          * Since this is only implemented for domain decomposition
992          * and domain decomposition does not use the graph,
993          * we do not need to worry about shifting.
994          */
995         gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
996                                  lambda[efptCOUL], lambda[efptVDW],
997                                  (stepWork.computeVirial || stepWork.computeEnergy),
998                                  step, wcycle);
999     }
1000 #endif /* GMX_MPI */
1001
1002     // Coordinates on the device are needed if PME or BufferOps are offloaded.
1003     // The local coordinates can be copied right away.
1004     // NOTE: Consider moving this copy to right after they are updated and constrained,
1005     //       if the later is not offloaded.
1006     if (useGpuPmeOnThisRank || useGpuXBufOps == BufferOpsUseGpu::True)
1007     {
1008         if (stepWork.doNeighborSearch)
1009         {
1010             stateGpu->reinit(mdatoms->homenr, cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1011             if (useGpuPmeOnThisRank)
1012             {
1013                 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1014                 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1015             }
1016         }
1017         stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), gmx::StatePropagatorDataGpu::AtomLocality::Local);
1018     }
1019
1020     if (useGpuPmeOnThisRank)
1021     {
1022         launchPmeGpuSpread(fr->pmedata, box, stepWork, pmeFlags, useGpuPmeFReduction, wcycle);
1023     }
1024
1025     /* do gridding for pair search */
1026     if (stepWork.doNeighborSearch)
1027     {
1028         if (graph && stepWork.stateChanged)
1029         {
1030             /* Calculate intramolecular shift vectors to make molecules whole */
1031             mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
1032         }
1033
1034         // TODO
1035         // - vzero is constant, do we need to pass it?
1036         // - box_diag should be passed directly to nbnxn_put_on_grid
1037         //
1038         rvec vzero;
1039         clear_rvec(vzero);
1040
1041         rvec box_diag;
1042         box_diag[XX] = box[XX][XX];
1043         box_diag[YY] = box[YY][YY];
1044         box_diag[ZZ] = box[ZZ][ZZ];
1045
1046         wallcycle_start(wcycle, ewcNS);
1047         if (!DOMAINDECOMP(cr))
1048         {
1049             wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1050             nbnxn_put_on_grid(nbv, box,
1051                               0, vzero, box_diag,
1052                               nullptr, 0, mdatoms->homenr, -1,
1053                               fr->cginfo, x.unpaddedArrayRef(),
1054                               0, nullptr);
1055             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1056         }
1057         else
1058         {
1059             wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1060             nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd),
1061                                        fr->cginfo, x.unpaddedArrayRef());
1062             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1063         }
1064
1065         nbv->setAtomProperties(*mdatoms, fr->cginfo);
1066
1067         wallcycle_stop(wcycle, ewcNS);
1068
1069         /* initialize the GPU nbnxm atom data and bonded data structures */
1070         if (simulationWork.useGpuNonbonded)
1071         {
1072             wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1073
1074             wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1075             Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1076             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1077
1078             if (fr->gpuBonded)
1079             {
1080                 /* Now we put all atoms on the grid, we can assign bonded
1081                  * interactions to the GPU, where the grid order is
1082                  * needed. Also the xq, f and fshift device buffers have
1083                  * been reallocated if needed, so the bonded code can
1084                  * learn about them. */
1085                 // TODO the xq, f, and fshift buffers are now shared
1086                 // resources, so they should be maintained by a
1087                 // higher-level object than the nb module.
1088                 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1089                                                                       top->idef,
1090                                                                       Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1091                                                                       Nbnxm::gpu_get_f(nbv->gpu_nbv),
1092                                                                       Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1093             }
1094             wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1095         }
1096     }
1097
1098     if (stepWork.doNeighborSearch)
1099     {
1100         // Need to run after the GPU-offload bonded interaction lists
1101         // are set up to be able to determine whether there is bonded work.
1102         runScheduleWork->domainWork =
1103             setupDomainLifetimeWorkload(*inputrec,
1104                                         *fr,
1105                                         pull_work,
1106                                         ed,
1107                                         top->idef,
1108                                         *fcd,
1109                                         *mdatoms,
1110                                         stepWork);
1111
1112         wallcycle_start_nocount(wcycle, ewcNS);
1113         wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1114         /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1115         nbv->constructPairlist(Nbnxm::InteractionLocality::Local,
1116                                &top->excls, step, nrnb);
1117
1118         nbv->setupGpuShortRangeWork(fr->gpuBonded, Nbnxm::InteractionLocality::Local);
1119
1120         wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1121         wallcycle_stop(wcycle, ewcNS);
1122
1123         if (useGpuXBufOps == BufferOpsUseGpu::True)
1124         {
1125             nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1126         }
1127         // For force buffer ops, we use the below conditon rather than
1128         // useGpuFBufOps to ensure that init is performed even if this
1129         // NS step is also a virial step (on which f buf ops are deactivated).
1130         if (simulationWork.useGpuBufferOps && simulationWork.useGpuNonbonded && (GMX_GPU == GMX_GPU_CUDA))
1131         {
1132             nbv->atomdata_init_add_nbat_f_to_f_gpu();
1133         }
1134     }
1135     else if (!EI_TPI(inputrec->eI))
1136     {
1137         if (useGpuXBufOps == BufferOpsUseGpu::True)
1138         {
1139             nbv->convertCoordinatesGpu(Nbnxm::AtomLocality::Local, false,
1140                                        stateGpu->getCoordinates());
1141         }
1142         else
1143         {
1144             nbv->convertCoordinates(Nbnxm::AtomLocality::Local, false,
1145                                     x.unpaddedArrayRef());
1146         }
1147     }
1148
1149     const gmx::DomainLifetimeWorkload &domainWork = runScheduleWork->domainWork;
1150
1151     if (simulationWork.useGpuNonbonded)
1152     {
1153         ddBalanceRegionHandler.openBeforeForceComputationGpu();
1154
1155         wallcycle_start(wcycle, ewcLAUNCH_GPU);
1156
1157         wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1158         Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1159         if (stepWork.doNeighborSearch || (useGpuXBufOps == BufferOpsUseGpu::False))
1160         {
1161             Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1162                                       Nbnxm::AtomLocality::Local);
1163         }
1164         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1165         // with X buffer ops offloaded to the GPU on all but the search steps
1166
1167         // bonded work not split into separate local and non-local, so with DD
1168         // we can only launch the kernel after non-local coordinates have been received.
1169         if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1170         {
1171             wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1172             fr->gpuBonded->launchKernel(fr, stepWork, box);
1173             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1174         }
1175
1176         /* launch local nonbonded work on GPU */
1177         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1178         do_nb_verlet(fr, ic, enerd, stepWork, Nbnxm::InteractionLocality::Local, enbvClearFNo,
1179                      step, nrnb, wcycle);
1180         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1181         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1182     }
1183
1184     if (useGpuPmeOnThisRank)
1185     {
1186         // In PME GPU and mixed mode we launch FFT / gather after the
1187         // X copy/transform to allow overlap as well as after the GPU NB
1188         // launch to avoid FFT launch overhead hijacking the CPU and delaying
1189         // the nonbonded kernel.
1190         launchPmeGpuFftAndGather(fr->pmedata, wcycle);
1191     }
1192
1193     // TODO Update this comment when introducing SimulationWorkload
1194     //
1195     // The conditions for gpuHaloExchange e.g. using GPU buffer
1196     // operations were checked before construction, so here we can
1197     // just use it and assert upon any conditions.
1198     gmx::GpuHaloExchange *gpuHaloExchange              = (havePPDomainDecomposition(cr) ? cr->dd->gpuHaloExchange.get() : nullptr);
1199     const bool            ddUsesGpuDirectCommunication = (gpuHaloExchange != nullptr);
1200     GMX_ASSERT(!ddUsesGpuDirectCommunication || (useGpuXBufOps == BufferOpsUseGpu::True),
1201                "Must use coordinate buffer ops with GPU halo exchange");
1202
1203     /* Communicate coordinates and sum dipole if necessary +
1204        do non-local pair search */
1205     if (havePPDomainDecomposition(cr))
1206     {
1207         if (stepWork.doNeighborSearch)
1208         {
1209             // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1210             wallcycle_start_nocount(wcycle, ewcNS);
1211             wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1212             /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1213             nbv->constructPairlist(Nbnxm::InteractionLocality::NonLocal,
1214                                    &top->excls, step, nrnb);
1215
1216             nbv->setupGpuShortRangeWork(fr->gpuBonded, Nbnxm::InteractionLocality::NonLocal);
1217             wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1218             wallcycle_stop(wcycle, ewcNS);
1219             if (ddUsesGpuDirectCommunication)
1220             {
1221                 gpuHaloExchange->reinitHalo(stateGpu->getCoordinates(), stateGpu->getForces());
1222             }
1223         }
1224         else
1225         {
1226             if (ddUsesGpuDirectCommunication)
1227             {
1228                 // The following must be called after local setCoordinates (which records an event
1229                 // when the coordinate data has been copied to the device).
1230                 gpuHaloExchange->communicateHaloCoordinates(box);
1231
1232                 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1233                 {
1234                     //non-local part of coordinate buffer must be copied back to host for CPU work
1235                     stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), gmx::StatePropagatorDataGpu::AtomLocality::NonLocal);
1236                 }
1237             }
1238             else
1239             {
1240                 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1241             }
1242
1243             if (useGpuXBufOps == BufferOpsUseGpu::True)
1244             {
1245                 // The condition here was (pme != nullptr && pme_gpu_get_device_x(fr->pmedata) != nullptr)
1246                 if (!useGpuPmeOnThisRank && !ddUsesGpuDirectCommunication)
1247                 {
1248                     stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), gmx::StatePropagatorDataGpu::AtomLocality::NonLocal);
1249                 }
1250                 nbv->convertCoordinatesGpu(Nbnxm::AtomLocality::NonLocal, false,
1251                                            stateGpu->getCoordinates());
1252             }
1253             else
1254             {
1255                 nbv->convertCoordinates(Nbnxm::AtomLocality::NonLocal, false,
1256                                         x.unpaddedArrayRef());
1257             }
1258
1259         }
1260
1261         if (simulationWork.useGpuNonbonded)
1262         {
1263             wallcycle_start(wcycle, ewcLAUNCH_GPU);
1264
1265             if (stepWork.doNeighborSearch || (useGpuXBufOps == BufferOpsUseGpu::False))
1266             {
1267                 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1268                 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1269                                           Nbnxm::AtomLocality::NonLocal);
1270                 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1271             }
1272
1273             if (domainWork.haveGpuBondedWork)
1274             {
1275                 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1276                 fr->gpuBonded->launchKernel(fr, stepWork, box);
1277                 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1278             }
1279
1280             /* launch non-local nonbonded tasks on GPU */
1281             wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1282             do_nb_verlet(fr, ic, enerd, stepWork, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1283                          step, nrnb, wcycle);
1284             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1285
1286             wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1287         }
1288     }
1289
1290     if (simulationWork.useGpuNonbonded)
1291     {
1292         /* launch D2H copy-back F */
1293         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1294         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1295
1296         bool copyBackNbForce  = (useGpuFBufOps == BufferOpsUseGpu::False);
1297
1298         if (havePPDomainDecomposition(cr))
1299         {
1300             Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1301                                       stepWork, Nbnxm::AtomLocality::NonLocal, copyBackNbForce);
1302         }
1303         Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1304                                   stepWork, Nbnxm::AtomLocality::Local, copyBackNbForce);
1305         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1306
1307         if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1308         {
1309             fr->gpuBonded->launchEnergyTransfer();
1310         }
1311         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1312     }
1313
1314     if (stepWork.stateChanged && inputrecNeedMutot(inputrec))
1315     {
1316         if (PAR(cr))
1317         {
1318             gmx_sumd(2*DIM, mu, cr);
1319
1320             ddBalanceRegionHandler.reopenRegionCpu();
1321         }
1322
1323         for (i = 0; i < 2; i++)
1324         {
1325             for (j = 0; j < DIM; j++)
1326             {
1327                 fr->mu_tot[i][j] = mu[i*DIM + j];
1328             }
1329         }
1330     }
1331     if (mdatoms->nChargePerturbed == 0)
1332     {
1333         copy_rvec(fr->mu_tot[0], mu_tot);
1334     }
1335     else
1336     {
1337         for (j = 0; j < DIM; j++)
1338         {
1339             mu_tot[j] =
1340                 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1341                 lambda[efptCOUL]*fr->mu_tot[1][j];
1342         }
1343     }
1344
1345     /* Reset energies */
1346     reset_enerdata(enerd);
1347     /* Clear the shift forces */
1348     // TODO: This should be linked to the shift force buffer in use, or cleared before use instead
1349     for (gmx::RVec &elem : fr->shiftForces)
1350     {
1351         elem = { 0.0_real, 0.0_real, 0.0_real };
1352     }
1353
1354     if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1355     {
1356         wallcycle_start(wcycle, ewcPPDURINGPME);
1357         dd_force_flop_start(cr->dd, nrnb);
1358     }
1359
1360     if (inputrec->bRot)
1361     {
1362         wallcycle_start(wcycle, ewcROT);
1363         do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, stepWork.doNeighborSearch);
1364         wallcycle_stop(wcycle, ewcROT);
1365     }
1366
1367     /* Start the force cycle counter.
1368      * Note that a different counter is used for dynamic load balancing.
1369      */
1370     wallcycle_start(wcycle, ewcFORCE);
1371
1372     // Set up and clear force outputs.
1373     // We use std::move to keep the compiler happy, it has no effect.
1374     ForceOutputs forceOut = setupForceOutputs(fr, pull_work, *inputrec, std::move(force), stepWork, wcycle);
1375
1376     /* We calculate the non-bonded forces, when done on the CPU, here.
1377      * We do this before calling do_force_lowlevel, because in that
1378      * function, the listed forces are calculated before PME, which
1379      * does communication.  With this order, non-bonded and listed
1380      * force calculation imbalance can be balanced out by the domain
1381      * decomposition load balancing.
1382      */
1383
1384     const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1385
1386     if (!useOrEmulateGpuNb)
1387     {
1388         do_nb_verlet(fr, ic, enerd, stepWork, Nbnxm::InteractionLocality::Local, enbvClearFYes,
1389                      step, nrnb, wcycle);
1390     }
1391
1392     if (fr->efep != efepNO)
1393     {
1394         /* Calculate the local and non-local free energy interactions here.
1395          * Happens here on the CPU both with and without GPU.
1396          */
1397         nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::Local,
1398                                       fr, as_rvec_array(x.unpaddedArrayRef().data()), &forceOut.forceWithShiftForces(), *mdatoms,
1399                                       inputrec->fepvals, lambda.data(),
1400                                       enerd, stepWork, nrnb);
1401
1402         if (havePPDomainDecomposition(cr))
1403         {
1404             nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::NonLocal,
1405                                           fr, as_rvec_array(x.unpaddedArrayRef().data()), &forceOut.forceWithShiftForces(), *mdatoms,
1406                                           inputrec->fepvals, lambda.data(),
1407                                           enerd, stepWork, nrnb);
1408         }
1409     }
1410
1411     if (!useOrEmulateGpuNb)
1412     {
1413         if (havePPDomainDecomposition(cr))
1414         {
1415             do_nb_verlet(fr, ic, enerd, stepWork, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1416                          step, nrnb, wcycle);
1417         }
1418
1419         if (stepWork.computeForces)
1420         {
1421             /* Add all the non-bonded force to the normal force array.
1422              * This can be split into a local and a non-local part when overlapping
1423              * communication with calculation with domain decomposition.
1424              */
1425             wallcycle_stop(wcycle, ewcFORCE);
1426             nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::All, forceOut.forceWithShiftForces().force());
1427             wallcycle_start_nocount(wcycle, ewcFORCE);
1428         }
1429
1430         /* If there are multiple fshift output buffers we need to reduce them */
1431         if (stepWork.computeVirial)
1432         {
1433             /* This is not in a subcounter because it takes a
1434                negligible and constant-sized amount of time */
1435             nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat,
1436                                                      forceOut.forceWithShiftForces().shiftForces());
1437         }
1438     }
1439
1440     /* update QMMMrec, if necessary */
1441     if (fr->bQMMM)
1442     {
1443         update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1444     }
1445
1446     // TODO Force flags should include haveFreeEnergyWork for this domain
1447     if (ddUsesGpuDirectCommunication &&
1448         (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1449     {
1450         /* Wait for non-local coordinate data to be copied from device */
1451         nbv->wait_nonlocal_x_copy_D2H_done();
1452     }
1453     /* Compute the bonded and non-bonded energies and optionally forces */
1454     do_force_lowlevel(fr, inputrec, &(top->idef),
1455                       cr, ms, nrnb, wcycle, mdatoms,
1456                       x, hist, &forceOut, enerd, fcd,
1457                       box, lambda.data(), graph, fr->mu_tot,
1458                       stepWork,
1459                       ddBalanceRegionHandler);
1460
1461     wallcycle_stop(wcycle, ewcFORCE);
1462
1463     computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1464                          imdSession, pull_work, step, t, wcycle,
1465                          fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda.data(),
1466                          stepWork, &forceOut.forceWithVirial(), enerd,
1467                          ed, stepWork.doNeighborSearch);
1468
1469
1470     // Will store the amount of cycles spent waiting for the GPU that
1471     // will be later used in the DLB accounting.
1472     float cycles_wait_gpu = 0;
1473     if (useOrEmulateGpuNb)
1474     {
1475         auto &forceWithShiftForces = forceOut.forceWithShiftForces();
1476
1477         /* wait for non-local forces (or calculate in emulation mode) */
1478         if (havePPDomainDecomposition(cr))
1479         {
1480             if (simulationWork.useGpuNonbonded)
1481             {
1482                 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1483                                                                stepWork, Nbnxm::AtomLocality::NonLocal,
1484                                                                enerd->grpp.ener[egLJSR].data(),
1485                                                                enerd->grpp.ener[egCOULSR].data(),
1486                                                                forceWithShiftForces.shiftForces(),
1487                                                                wcycle);
1488             }
1489             else
1490             {
1491                 wallcycle_start_nocount(wcycle, ewcFORCE);
1492                 do_nb_verlet(fr, ic, enerd, stepWork, Nbnxm::InteractionLocality::NonLocal, enbvClearFYes,
1493                              step, nrnb, wcycle);
1494                 wallcycle_stop(wcycle, ewcFORCE);
1495             }
1496
1497             if (useGpuFBufOps == BufferOpsUseGpu::True)
1498             {
1499                 gmx::FixedCapacityVector<GpuEventSynchronizer*, 1> dependencyList;
1500
1501                 // TODO: move this into DomainLifetimeWorkload, including the second part of the condition
1502                 // The bonded and free energy CPU tasks can have non-local force contributions
1503                 // which are a dependency for the GPU force reduction.
1504                 bool  haveNonLocalForceContribInCpuBuffer = domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1505
1506                 if (haveNonLocalForceContribInCpuBuffer)
1507                 {
1508                     stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(), gmx::StatePropagatorDataGpu::AtomLocality::NonLocal);
1509                     dependencyList.push_back(stateGpu->getForcesReadyOnDeviceEvent(gmx::StatePropagatorDataGpu::AtomLocality::NonLocal,
1510                                                                                    useGpuFBufOps == BufferOpsUseGpu::True));
1511                 }
1512
1513                 nbv->atomdata_add_nbat_f_to_f_gpu(Nbnxm::AtomLocality::NonLocal,
1514                                                   stateGpu->getForces(),
1515                                                   pme_gpu_get_device_f(fr->pmedata),
1516                                                   dependencyList,
1517                                                   false, haveNonLocalForceContribInCpuBuffer);
1518                 stateGpu->copyForcesFromGpu(forceOut.forceWithShiftForces().force(), gmx::StatePropagatorDataGpu::AtomLocality::NonLocal);
1519             }
1520             else
1521             {
1522                 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::NonLocal,
1523                                               forceWithShiftForces.force());
1524             }
1525
1526
1527             if (fr->nbv->emulateGpu() && stepWork.computeVirial)
1528             {
1529                 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat,
1530                                                          forceWithShiftForces.shiftForces());
1531             }
1532         }
1533     }
1534
1535     const bool useGpuForcesHaloExchange = ddUsesGpuDirectCommunication && (useGpuFBufOps == BufferOpsUseGpu::True);
1536     const bool useCpuPmeFReduction      = thisRankHasDuty(cr, DUTY_PME) && !useGpuPmeFReduction;
1537     // TODO: move this into DomainLifetimeWorkload, including the second part of the condition
1538     const bool haveCpuLocalForces     = (domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork || useCpuPmeFReduction ||
1539                                          (fr->efep != efepNO));
1540
1541     if (havePPDomainDecomposition(cr))
1542     {
1543         /* We are done with the CPU compute.
1544          * We will now communicate the non-local forces.
1545          * If we use a GPU this will overlap with GPU work, so in that case
1546          * we do not close the DD force balancing region here.
1547          */
1548         ddBalanceRegionHandler.closeAfterForceComputationCpu();
1549
1550         if (stepWork.computeForces)
1551         {
1552
1553             if (useGpuForcesHaloExchange)
1554             {
1555                 if (haveCpuLocalForces)
1556                 {
1557                     stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(), gmx::StatePropagatorDataGpu::AtomLocality::Local);
1558                 }
1559                 gpuHaloExchange->communicateHaloForces(haveCpuLocalForces);
1560             }
1561             else
1562             {
1563                 if (useGpuFBufOps == BufferOpsUseGpu::True)
1564                 {
1565                     nbv->wait_for_gpu_force_reduction(Nbnxm::AtomLocality::NonLocal);
1566                 }
1567                 dd_move_f(cr->dd, &forceOut.forceWithShiftForces(), wcycle);
1568             }
1569
1570         }
1571     }
1572
1573     // With both nonbonded and PME offloaded a GPU on the same rank, we use
1574     // an alternating wait/reduction scheme.
1575     bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded && !DOMAINDECOMP(cr) &&
1576                              (useGpuFBufOps == BufferOpsUseGpu::False));
1577     if (alternateGpuWait)
1578     {
1579         alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &forceOut, enerd,
1580                                     stepWork, pmeFlags, wcycle);
1581     }
1582
1583     if (!alternateGpuWait && useGpuPmeOnThisRank)
1584     {
1585         pme_gpu_wait_and_reduce(fr->pmedata, pmeFlags, wcycle, &forceOut.forceWithVirial(), enerd);
1586     }
1587
1588     /* Wait for local GPU NB outputs on the non-alternating wait path */
1589     if (!alternateGpuWait && simulationWork.useGpuNonbonded)
1590     {
1591         /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1592          * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1593          * but even with a step of 0.1 ms the difference is less than 1%
1594          * of the step time.
1595          */
1596         const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
1597         const float waitCycles               =
1598             Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1599                                         stepWork, Nbnxm::AtomLocality::Local,
1600                                         enerd->grpp.ener[egLJSR].data(),
1601                                         enerd->grpp.ener[egCOULSR].data(),
1602                                         forceOut.forceWithShiftForces().shiftForces(),
1603                                         wcycle);
1604
1605         if (ddBalanceRegionHandler.useBalancingRegion())
1606         {
1607             DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1608             if (stepWork.computeForces &&  waitCycles <= gpuWaitApiOverheadMargin)
1609             {
1610                 /* We measured few cycles, it could be that the kernel
1611                  * and transfer finished earlier and there was no actual
1612                  * wait time, only API call overhead.
1613                  * Then the actual time could be anywhere between 0 and
1614                  * cycles_wait_est. We will use half of cycles_wait_est.
1615                  */
1616                 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1617             }
1618             ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1619         }
1620     }
1621
1622     if (fr->nbv->emulateGpu())
1623     {
1624         // NOTE: emulation kernel is not included in the balancing region,
1625         // but emulation mode does not target performance anyway
1626         wallcycle_start_nocount(wcycle, ewcFORCE);
1627         do_nb_verlet(fr, ic, enerd, stepWork, Nbnxm::InteractionLocality::Local,
1628                      DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1629                      step, nrnb, wcycle);
1630         wallcycle_stop(wcycle, ewcFORCE);
1631     }
1632
1633     /* Do the nonbonded GPU (or emulation) force buffer reduction
1634      * on the non-alternating path. */
1635     if (useOrEmulateGpuNb && !alternateGpuWait)
1636     {
1637         gmx::FixedCapacityVector<GpuEventSynchronizer*, 2> dependencyList;
1638
1639         if (useGpuPmeFReduction)
1640         {
1641             dependencyList.push_back(pme_gpu_get_f_ready_synchronizer(fr->pmedata));
1642         }
1643
1644         gmx::ArrayRef<gmx::RVec>  forceWithShift = forceOut.forceWithShiftForces().force();
1645
1646         if (useGpuFBufOps == BufferOpsUseGpu::True)
1647         {
1648             // Flag to specify whether the CPU force buffer has contributions to
1649             // local atoms. This depends on whether there are CPU-based force tasks
1650             // or when DD is active the halo exchange has resulted in contributions
1651             // from the non-local part.
1652             const bool haveLocalForceContribInCpuBuffer = (haveCpuLocalForces || havePPDomainDecomposition(cr));
1653
1654             // TODO: move these steps as early as possible:
1655             // - CPU f H2D should be as soon as all CPU-side forces are done
1656             // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
1657             //   before the next CPU task that consumes the forces: vsite spread or update)
1658             // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
1659             //   of the halo exchange. In that case the copy is instead performed above, before the exchange.
1660             //   These should be unified.
1661             if (haveLocalForceContribInCpuBuffer && !useGpuForcesHaloExchange)
1662             {
1663                 stateGpu->copyForcesToGpu(forceWithShift, gmx::StatePropagatorDataGpu::AtomLocality::Local);
1664                 dependencyList.push_back(stateGpu->getForcesReadyOnDeviceEvent(gmx::StatePropagatorDataGpu::AtomLocality::Local,
1665                                                                                useGpuFBufOps == BufferOpsUseGpu::True));
1666             }
1667             if (useGpuForcesHaloExchange)
1668             {
1669                 // Add a stream synchronization to satisfy a dependency
1670                 // for the local buffer ops on the result of GPU halo
1671                 // exchange, which operates in the non-local stream and
1672                 // writes to to local parf og the force buffer.
1673                 //
1674                 // TODO improve this through use of an event - see Redmine #3093
1675                 //      push the event into the dependencyList
1676                 nbv->stream_local_wait_for_nonlocal();
1677             }
1678             nbv->atomdata_add_nbat_f_to_f_gpu(Nbnxm::AtomLocality::Local,
1679                                               stateGpu->getForces(),
1680                                               pme_gpu_get_device_f(fr->pmedata),
1681                                               dependencyList,
1682                                               useGpuPmeFReduction, haveLocalForceContribInCpuBuffer);
1683             // This function call synchronizes the local stream
1684             nbv->wait_for_gpu_force_reduction(Nbnxm::AtomLocality::Local);
1685             stateGpu->copyForcesFromGpu(forceWithShift, gmx::StatePropagatorDataGpu::AtomLocality::Local);
1686         }
1687         else
1688         {
1689             nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local, forceWithShift);
1690         }
1691
1692     }
1693
1694     launchGpuEndOfStepTasks(nbv, fr->gpuBonded, fr->pmedata, enerd,
1695                             *runScheduleWork,
1696                             simulationWork.useGpuNonbonded, useGpuPmeOnThisRank,
1697                             step,
1698                             wcycle);
1699
1700     if (DOMAINDECOMP(cr))
1701     {
1702         dd_force_flop_stop(cr->dd, nrnb);
1703     }
1704
1705     if (stepWork.computeForces)
1706     {
1707         rvec *f = as_rvec_array(forceOut.forceWithShiftForces().force().data());
1708
1709         /* If we have NoVirSum forces, but we do not calculate the virial,
1710          * we sum fr->f_novirsum=forceOut.f later.
1711          */
1712         if (vsite && !(fr->haveDirectVirialContributions && !stepWork.computeVirial))
1713         {
1714             rvec *fshift = as_rvec_array(forceOut.forceWithShiftForces().shiftForces().data());
1715             spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), f, fshift, FALSE, nullptr, nrnb,
1716                            &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1717         }
1718
1719         if (stepWork.computeVirial)
1720         {
1721             /* Calculation of the virial must be done after vsites! */
1722             calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()),
1723                         forceOut.forceWithShiftForces(),
1724                         vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1725         }
1726     }
1727
1728     if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1729     {
1730         /* In case of node-splitting, the PP nodes receive the long-range
1731          * forces, virial and energy from the PME nodes here.
1732          */
1733         pme_receive_force_ener(cr, &forceOut.forceWithVirial(), enerd, wcycle);
1734     }
1735
1736     if (stepWork.computeForces)
1737     {
1738         post_process_forces(cr, step, nrnb, wcycle,
1739                             top, box, as_rvec_array(x.unpaddedArrayRef().data()), &forceOut,
1740                             vir_force, mdatoms, graph, fr, vsite,
1741                             stepWork);
1742     }
1743
1744     if (stepWork.computeEnergy)
1745     {
1746         /* Sum the potential energy terms from group contributions */
1747         sum_epot(&(enerd->grpp), enerd->term);
1748
1749         if (!EI_TPI(inputrec->eI))
1750         {
1751             checkPotentialEnergyValidity(step, *enerd, *inputrec);
1752         }
1753     }
1754
1755     /* In case we don't have constraints and are using GPUs, the next balancing
1756      * region starts here.
1757      * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1758      * virial calculation and COM pulling, is not thus not included in
1759      * the balance timing, which is ok as most tasks do communication.
1760      */
1761     ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);
1762 }