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