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