Use ArrayRef in signatures of constraints and some pulling
[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,2021, 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 #include <optional>
48
49 #include "gromacs/applied_forces/awh/awh.h"
50 #include "gromacs/domdec/dlbtiming.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/domdec/gpuhaloexchange.h"
54 #include "gromacs/domdec/partition.h"
55 #include "gromacs/essentialdynamics/edsam.h"
56 #include "gromacs/ewald/pme.h"
57 #include "gromacs/ewald/pme_pp.h"
58 #include "gromacs/ewald/pme_pp_comm_gpu.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
61 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
62 #include "gromacs/gmxlib/nrnb.h"
63 #include "gromacs/gpu_utils/gpu_utils.h"
64 #include "gromacs/imd/imd.h"
65 #include "gromacs/listed_forces/disre.h"
66 #include "gromacs/listed_forces/gpubonded.h"
67 #include "gromacs/listed_forces/listed_forces.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/dispersioncorrection.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/force.h"
80 #include "gromacs/mdlib/force_flags.h"
81 #include "gromacs/mdlib/forcerec.h"
82 #include "gromacs/mdlib/gmx_omp_nthreads.h"
83 #include "gromacs/mdlib/update.h"
84 #include "gromacs/mdlib/vsite.h"
85 #include "gromacs/mdlib/wall.h"
86 #include "gromacs/mdlib/wholemoleculetransform.h"
87 #include "gromacs/mdtypes/commrec.h"
88 #include "gromacs/mdtypes/enerdata.h"
89 #include "gromacs/mdtypes/forcebuffers.h"
90 #include "gromacs/mdtypes/forceoutput.h"
91 #include "gromacs/mdtypes/forcerec.h"
92 #include "gromacs/mdtypes/iforceprovider.h"
93 #include "gromacs/mdtypes/inputrec.h"
94 #include "gromacs/mdtypes/md_enums.h"
95 #include "gromacs/mdtypes/mdatom.h"
96 #include "gromacs/mdtypes/multipletimestepping.h"
97 #include "gromacs/mdtypes/simulation_workload.h"
98 #include "gromacs/mdtypes/state.h"
99 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
100 #include "gromacs/nbnxm/gpu_data_mgmt.h"
101 #include "gromacs/nbnxm/nbnxm.h"
102 #include "gromacs/nbnxm/nbnxm_gpu.h"
103 #include "gromacs/pbcutil/ishift.h"
104 #include "gromacs/pbcutil/pbc.h"
105 #include "gromacs/pulling/pull.h"
106 #include "gromacs/pulling/pull_rotation.h"
107 #include "gromacs/timing/cyclecounter.h"
108 #include "gromacs/timing/gpu_timing.h"
109 #include "gromacs/timing/wallcycle.h"
110 #include "gromacs/timing/wallcyclereporting.h"
111 #include "gromacs/timing/walltime_accounting.h"
112 #include "gromacs/topology/topology.h"
113 #include "gromacs/utility/arrayref.h"
114 #include "gromacs/utility/basedefinitions.h"
115 #include "gromacs/utility/cstringutil.h"
116 #include "gromacs/utility/exceptions.h"
117 #include "gromacs/utility/fatalerror.h"
118 #include "gromacs/utility/fixedcapacityvector.h"
119 #include "gromacs/utility/gmxassert.h"
120 #include "gromacs/utility/gmxmpi.h"
121 #include "gromacs/utility/logger.h"
122 #include "gromacs/utility/smalloc.h"
123 #include "gromacs/utility/strconvert.h"
124 #include "gromacs/utility/sysinfo.h"
125
126 #include "gpuforcereduction.h"
127
128 using gmx::ArrayRef;
129 using gmx::AtomLocality;
130 using gmx::DomainLifetimeWorkload;
131 using gmx::ForceOutputs;
132 using gmx::ForceWithShiftForces;
133 using gmx::InteractionLocality;
134 using gmx::RVec;
135 using gmx::SimulationWorkload;
136 using gmx::StepWorkload;
137
138 // TODO: this environment variable allows us to verify before release
139 // that on less common architectures the total cost of polling is not larger than
140 // a blocking wait (so polling does not introduce overhead when the static
141 // PME-first ordering would suffice).
142 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
143
144 static void sum_forces(ArrayRef<RVec> f, ArrayRef<const RVec> forceToAdd)
145 {
146     GMX_ASSERT(f.size() >= forceToAdd.size(), "Accumulation buffer should be sufficiently large");
147     const int end = forceToAdd.size();
148
149     int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
150 #pragma omp parallel for num_threads(nt) schedule(static)
151     for (int i = 0; i < end; i++)
152     {
153         rvec_inc(f[i], forceToAdd[i]);
154     }
155 }
156
157 static void calc_virial(int                              start,
158                         int                              homenr,
159                         const rvec                       x[],
160                         const gmx::ForceWithShiftForces& forceWithShiftForces,
161                         tensor                           vir_part,
162                         const matrix                     box,
163                         t_nrnb*                          nrnb,
164                         const t_forcerec*                fr,
165                         PbcType                          pbcType)
166 {
167     /* The short-range virial from surrounding boxes */
168     const rvec* fshift = as_rvec_array(forceWithShiftForces.shiftForces().data());
169     calc_vir(SHIFTS, fr->shift_vec, fshift, vir_part, pbcType == PbcType::Screw, box);
170     inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
171
172     /* Calculate partial virial, for local atoms only, based on short range.
173      * Total virial is computed in global_stat, called from do_md
174      */
175     const rvec* f = as_rvec_array(forceWithShiftForces.force().data());
176     f_calc_vir(start, start + homenr, x, f, vir_part, box);
177     inc_nrnb(nrnb, eNR_VIRIAL, homenr);
178
179     if (debug)
180     {
181         pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
182     }
183 }
184
185 static void pull_potential_wrapper(const t_commrec*               cr,
186                                    const t_inputrec&              ir,
187                                    const matrix                   box,
188                                    gmx::ArrayRef<const gmx::RVec> x,
189                                    gmx::ForceWithVirial*          force,
190                                    const t_mdatoms*               mdatoms,
191                                    gmx_enerdata_t*                enerd,
192                                    pull_t*                        pull_work,
193                                    const real*                    lambda,
194                                    double                         t,
195                                    gmx_wallcycle_t                wcycle)
196 {
197     t_pbc pbc;
198     real  dvdl;
199
200     /* Calculate the center of mass forces, this requires communication,
201      * which is why pull_potential is called close to other communication.
202      */
203     wallcycle_start(wcycle, ewcPULLPOT);
204     set_pbc(&pbc, ir.pbcType, box);
205     dvdl = 0;
206     enerd->term[F_COM_PULL] +=
207             pull_potential(pull_work,
208                            gmx::arrayRefFromArray(mdatoms->massT, mdatoms->nr),
209                            &pbc,
210                            cr,
211                            t,
212                            lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Restraint)],
213                            x,
214                            force,
215                            &dvdl);
216     enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Restraint] += dvdl;
217     wallcycle_stop(wcycle, ewcPULLPOT);
218 }
219
220 static void pme_receive_force_ener(t_forcerec*           fr,
221                                    const t_commrec*      cr,
222                                    gmx::ForceWithVirial* forceWithVirial,
223                                    gmx_enerdata_t*       enerd,
224                                    bool                  useGpuPmePpComms,
225                                    bool                  receivePmeForceToGpu,
226                                    gmx_wallcycle_t       wcycle)
227 {
228     real  e_q, e_lj, dvdl_q, dvdl_lj;
229     float cycles_ppdpme, cycles_seppme;
230
231     cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
232     dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
233
234     /* In case of node-splitting, the PP nodes receive the long-range
235      * forces, virial and energy from the PME nodes here.
236      */
237     wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
238     dvdl_q  = 0;
239     dvdl_lj = 0;
240     gmx_pme_receive_f(fr->pmePpCommGpu.get(),
241                       cr,
242                       forceWithVirial,
243                       &e_q,
244                       &e_lj,
245                       &dvdl_q,
246                       &dvdl_lj,
247                       useGpuPmePpComms,
248                       receivePmeForceToGpu,
249                       &cycles_seppme);
250     enerd->term[F_COUL_RECIP] += e_q;
251     enerd->term[F_LJ_RECIP] += e_lj;
252     enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Coul] += dvdl_q;
253     enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += dvdl_lj;
254
255     if (wcycle)
256     {
257         dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
258     }
259     wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
260 }
261
262 static void print_large_forces(FILE*                fp,
263                                const t_mdatoms*     md,
264                                const t_commrec*     cr,
265                                int64_t              step,
266                                real                 forceTolerance,
267                                ArrayRef<const RVec> x,
268                                ArrayRef<const RVec> f)
269 {
270     real       force2Tolerance = gmx::square(forceTolerance);
271     gmx::index numNonFinite    = 0;
272     for (int i = 0; i < md->homenr; i++)
273     {
274         real force2    = norm2(f[i]);
275         bool nonFinite = !std::isfinite(force2);
276         if (force2 >= force2Tolerance || nonFinite)
277         {
278             fprintf(fp,
279                     "step %" PRId64 " atom %6d  x %8.3f %8.3f %8.3f  force %12.5e\n",
280                     step,
281                     ddglatnr(cr->dd, i),
282                     x[i][XX],
283                     x[i][YY],
284                     x[i][ZZ],
285                     std::sqrt(force2));
286         }
287         if (nonFinite)
288         {
289             numNonFinite++;
290         }
291     }
292     if (numNonFinite > 0)
293     {
294         /* Note that with MPI this fatal call on one rank might interrupt
295          * the printing on other ranks. But we can only avoid that with
296          * an expensive MPI barrier that we would need at each step.
297          */
298         gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
299     }
300 }
301
302 //! When necessary, spreads forces on vsites and computes the virial for \p forceOutputs->forceWithShiftForces()
303 static void postProcessForceWithShiftForces(t_nrnb*                   nrnb,
304                                             gmx_wallcycle_t           wcycle,
305                                             const matrix              box,
306                                             ArrayRef<const RVec>      x,
307                                             ForceOutputs*             forceOutputs,
308                                             tensor                    vir_force,
309                                             const t_mdatoms&          mdatoms,
310                                             const t_forcerec&         fr,
311                                             gmx::VirtualSitesHandler* vsite,
312                                             const StepWorkload&       stepWork)
313 {
314     ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
315
316     /* If we have NoVirSum forces, but we do not calculate the virial,
317      * we later sum the forceWithShiftForces buffer together with
318      * the noVirSum buffer and spread the combined vsite forces at once.
319      */
320     if (vsite && (!forceOutputs->haveForceWithVirial() || stepWork.computeVirial))
321     {
322         using VirialHandling = gmx::VirtualSitesHandler::VirialHandling;
323
324         auto                 f      = forceWithShiftForces.force();
325         auto                 fshift = forceWithShiftForces.shiftForces();
326         const VirialHandling virialHandling =
327                 (stepWork.computeVirial ? VirialHandling::Pbc : VirialHandling::None);
328         vsite->spreadForces(x, f, virialHandling, fshift, nullptr, nrnb, box, wcycle);
329         forceWithShiftForces.haveSpreadVsiteForces() = true;
330     }
331
332     if (stepWork.computeVirial)
333     {
334         /* Calculation of the virial must be done after vsites! */
335         calc_virial(
336                 0, mdatoms.homenr, as_rvec_array(x.data()), forceWithShiftForces, vir_force, box, nrnb, &fr, fr.pbcType);
337     }
338 }
339
340 //! Spread, compute virial for and sum forces, when necessary
341 static void postProcessForces(const t_commrec*          cr,
342                               int64_t                   step,
343                               t_nrnb*                   nrnb,
344                               gmx_wallcycle_t           wcycle,
345                               const matrix              box,
346                               ArrayRef<const RVec>      x,
347                               ForceOutputs*             forceOutputs,
348                               tensor                    vir_force,
349                               const t_mdatoms*          mdatoms,
350                               const t_forcerec*         fr,
351                               gmx::VirtualSitesHandler* vsite,
352                               const StepWorkload&       stepWork)
353 {
354     // Extract the final output force buffer, which is also the buffer for forces with shift forces
355     ArrayRef<RVec> f = forceOutputs->forceWithShiftForces().force();
356
357     if (forceOutputs->haveForceWithVirial())
358     {
359         auto& forceWithVirial = forceOutputs->forceWithVirial();
360
361         if (vsite)
362         {
363             /* Spread the mesh force on virtual sites to the other particles...
364              * This is parallellized. MPI communication is performed
365              * if the constructing atoms aren't local.
366              */
367             GMX_ASSERT(!stepWork.computeVirial || f.data() != forceWithVirial.force_.data(),
368                        "We need separate force buffers for shift and virial forces when "
369                        "computing the virial");
370             GMX_ASSERT(!stepWork.computeVirial
371                                || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
372                        "We should spread the force with shift forces separately when computing "
373                        "the virial");
374             const gmx::VirtualSitesHandler::VirialHandling virialHandling =
375                     (stepWork.computeVirial ? gmx::VirtualSitesHandler::VirialHandling::NonLinear
376                                             : gmx::VirtualSitesHandler::VirialHandling::None);
377             matrix virial = { { 0 } };
378             vsite->spreadForces(x, forceWithVirial.force_, virialHandling, {}, virial, nrnb, box, wcycle);
379             forceWithVirial.addVirialContribution(virial);
380         }
381
382         if (stepWork.computeVirial)
383         {
384             /* Now add the forces, this is local */
385             sum_forces(f, forceWithVirial.force_);
386
387             /* Add the direct virial contributions */
388             GMX_ASSERT(
389                     forceWithVirial.computeVirial_,
390                     "forceWithVirial should request virial computation when we request the virial");
391             m_add(vir_force, forceWithVirial.getVirial(), vir_force);
392
393             if (debug)
394             {
395                 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
396             }
397         }
398     }
399     else
400     {
401         GMX_ASSERT(vsite == nullptr || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
402                    "We should have spread the vsite forces (earlier)");
403     }
404
405     if (fr->print_force >= 0)
406     {
407         print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
408     }
409 }
410
411 static void do_nb_verlet(t_forcerec*                fr,
412                          const interaction_const_t* ic,
413                          gmx_enerdata_t*            enerd,
414                          const StepWorkload&        stepWork,
415                          const InteractionLocality  ilocality,
416                          const int                  clearF,
417                          const int64_t              step,
418                          t_nrnb*                    nrnb,
419                          gmx_wallcycle_t            wcycle)
420 {
421     if (!stepWork.computeNonbondedForces)
422     {
423         /* skip non-bonded calculation */
424         return;
425     }
426
427     nonbonded_verlet_t* nbv = fr->nbv.get();
428
429     /* GPU kernel launch overhead is already timed separately */
430     if (!nbv->useGpu())
431     {
432         /* When dynamic pair-list  pruning is requested, we need to prune
433          * at nstlistPrune steps.
434          */
435         if (nbv->isDynamicPruningStepCpu(step))
436         {
437             /* Prune the pair-list beyond fr->ic->rlistPrune using
438              * the current coordinates of the atoms.
439              */
440             wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
441             nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
442             wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
443         }
444     }
445
446     nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
447 }
448
449 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
450 {
451     int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, v.ssize());
452
453     /* Note that we would like to avoid this conditional by putting it
454      * into the omp pragma instead, but then we still take the full
455      * omp parallel for overhead (at least with gcc5).
456      */
457     if (!useOpenmpThreading || nth == 1)
458     {
459         for (RVec& elem : v)
460         {
461             clear_rvec(elem);
462         }
463     }
464     else
465     {
466 #pragma omp parallel for num_threads(nth) schedule(static)
467         for (gmx::index i = 0; i < v.ssize(); i++)
468         {
469             clear_rvec(v[i]);
470         }
471     }
472 }
473
474 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
475  *
476  * \param groupOptions  Group options, containing T-coupling options
477  */
478 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
479 {
480     real nrdfCoupled   = 0;
481     real nrdfUncoupled = 0;
482     real kineticEnergy = 0;
483     for (int g = 0; g < groupOptions.ngtc; g++)
484     {
485         if (groupOptions.tau_t[g] >= 0)
486         {
487             nrdfCoupled += groupOptions.nrdf[g];
488             kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * gmx::c_boltz;
489         }
490         else
491         {
492             nrdfUncoupled += groupOptions.nrdf[g];
493         }
494     }
495
496     /* This conditional with > also catches nrdf=0 */
497     if (nrdfCoupled > nrdfUncoupled)
498     {
499         return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
500     }
501     else
502     {
503         return 0;
504     }
505 }
506
507 /*! \brief This routine checks that the potential energy is finite.
508  *
509  * Always checks that the potential energy is finite. If step equals
510  * inputrec.init_step also checks that the magnitude of the potential energy
511  * is reasonable. Terminates with a fatal error when a check fails.
512  * Note that passing this check does not guarantee finite forces,
513  * since those use slightly different arithmetics. But in most cases
514  * there is just a narrow coordinate range where forces are not finite
515  * and energies are finite.
516  *
517  * \param[in] step      The step number, used for checking and printing
518  * \param[in] enerd     The energy data; the non-bonded group energies need to be added to
519  * enerd.term[F_EPOT] before calling this routine \param[in] inputrec  The input record
520  */
521 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
522 {
523     /* Threshold valid for comparing absolute potential energy against
524      * the kinetic energy. Normally one should not consider absolute
525      * potential energy values, but with a factor of one million
526      * we should never get false positives.
527      */
528     constexpr real c_thresholdFactor = 1e6;
529
530     bool energyIsNotFinite    = !std::isfinite(enerd.term[F_EPOT]);
531     real averageKineticEnergy = 0;
532     /* We only check for large potential energy at the initial step,
533      * because that is by far the most likely step for this too occur
534      * and because computing the average kinetic energy is not free.
535      * Note: nstcalcenergy >> 1 often does not allow to catch large energies
536      * before they become NaN.
537      */
538     if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
539     {
540         averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
541     }
542
543     if (energyIsNotFinite
544         || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
545     {
546         gmx_fatal(
547                 FARGS,
548                 "Step %" PRId64
549                 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
550                 "contributions to the energy are %g and %g, respectively. A %s potential energy "
551                 "can be caused by overlapping interactions in bonded interactions or very large%s "
552                 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
553                 "configuration, incorrect interactions or parameters in the topology.",
554                 step,
555                 enerd.term[F_EPOT],
556                 energyIsNotFinite ? "not finite" : "extremely high",
557                 enerd.term[F_LJ],
558                 enerd.term[F_COUL_SR],
559                 energyIsNotFinite ? "non-finite" : "very high",
560                 energyIsNotFinite ? " or Nan" : "");
561     }
562 }
563
564 /*! \brief Return true if there are special forces computed this step.
565  *
566  * The conditionals exactly correspond to those in computeSpecialForces().
567  */
568 static bool haveSpecialForces(const t_inputrec&          inputrec,
569                               const gmx::ForceProviders& forceProviders,
570                               const pull_t*              pull_work,
571                               const bool                 computeForces,
572                               const gmx_edsam*           ed)
573 {
574
575     return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
576             (inputrec.bPull && pull_have_potential(*pull_work)) ||  // pull
577             inputrec.bRot ||                                        // enforced rotation
578             (ed != nullptr) ||                                      // flooding
579             (inputrec.bIMD && computeForces));                      // IMD
580 }
581
582 /*! \brief Compute forces and/or energies for special algorithms
583  *
584  * The intention is to collect all calls to algorithms that compute
585  * forces on local atoms only and that do not contribute to the local
586  * virial sum (but add their virial contribution separately).
587  * Eventually these should likely all become ForceProviders.
588  * Within this function the intention is to have algorithms that do
589  * global communication at the end, so global barriers within the MD loop
590  * are as close together as possible.
591  *
592  * \param[in]     fplog            The log file
593  * \param[in]     cr               The communication record
594  * \param[in]     inputrec         The input record
595  * \param[in]     awh              The Awh module (nullptr if none in use).
596  * \param[in]     enforcedRotation Enforced rotation module.
597  * \param[in]     imdSession       The IMD session
598  * \param[in]     pull_work        The pull work structure.
599  * \param[in]     step             The current MD step
600  * \param[in]     t                The current time
601  * \param[in,out] wcycle           Wallcycle accounting struct
602  * \param[in,out] forceProviders   Pointer to a list of force providers
603  * \param[in]     box              The unit cell
604  * \param[in]     x                The coordinates
605  * \param[in]     mdatoms          Per atom properties
606  * \param[in]     lambda           Array of free-energy lambda values
607  * \param[in]     stepWork         Step schedule flags
608  * \param[in,out] forceWithVirialMtsLevel0  Force and virial for MTS level0 forces
609  * \param[in,out] forceWithVirialMtsLevel1  Force and virial for MTS level1 forces, can be nullptr
610  * \param[in,out] enerd            Energy buffer
611  * \param[in,out] ed               Essential dynamics pointer
612  * \param[in]     didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
613  *
614  * \todo Remove didNeighborSearch, which is used incorrectly.
615  * \todo Convert all other algorithms called here to ForceProviders.
616  */
617 static void computeSpecialForces(FILE*                          fplog,
618                                  const t_commrec*               cr,
619                                  const t_inputrec&              inputrec,
620                                  gmx::Awh*                      awh,
621                                  gmx_enfrot*                    enforcedRotation,
622                                  gmx::ImdSession*               imdSession,
623                                  pull_t*                        pull_work,
624                                  int64_t                        step,
625                                  double                         t,
626                                  gmx_wallcycle_t                wcycle,
627                                  gmx::ForceProviders*           forceProviders,
628                                  const matrix                   box,
629                                  gmx::ArrayRef<const gmx::RVec> x,
630                                  const t_mdatoms*               mdatoms,
631                                  gmx::ArrayRef<const real>      lambda,
632                                  const StepWorkload&            stepWork,
633                                  gmx::ForceWithVirial*          forceWithVirialMtsLevel0,
634                                  gmx::ForceWithVirial*          forceWithVirialMtsLevel1,
635                                  gmx_enerdata_t*                enerd,
636                                  gmx_edsam*                     ed,
637                                  bool                           didNeighborSearch)
638 {
639     /* NOTE: Currently all ForceProviders only provide forces.
640      *       When they also provide energies, remove this conditional.
641      */
642     if (stepWork.computeForces)
643     {
644         gmx::ForceProviderInput forceProviderInput(
645                 x,
646                 mdatoms->homenr,
647                 gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->homenr),
648                 gmx::arrayRefFromArray(mdatoms->massT, mdatoms->homenr),
649                 t,
650                 box,
651                 *cr);
652         gmx::ForceProviderOutput forceProviderOutput(forceWithVirialMtsLevel0, enerd);
653
654         /* Collect forces from modules */
655         forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
656     }
657
658     if (inputrec.bPull && pull_have_potential(*pull_work))
659     {
660         const int mtsLevel = forceGroupMtsLevel(inputrec.mtsLevels, gmx::MtsForceGroups::Pull);
661         if (mtsLevel == 0 || stepWork.computeSlowForces)
662         {
663             auto& forceWithVirial = (mtsLevel == 0) ? forceWithVirialMtsLevel0 : forceWithVirialMtsLevel1;
664             pull_potential_wrapper(
665                     cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work, lambda.data(), t, wcycle);
666         }
667     }
668     if (awh)
669     {
670         const int mtsLevel = forceGroupMtsLevel(inputrec.mtsLevels, gmx::MtsForceGroups::Pull);
671         if (mtsLevel == 0 || stepWork.computeSlowForces)
672         {
673             const bool needForeignEnergyDifferences = awh->needForeignEnergyDifferences(step);
674             std::vector<double> foreignLambdaDeltaH, foreignLambdaDhDl;
675             if (needForeignEnergyDifferences)
676             {
677                 enerd->foreignLambdaTerms.finalizePotentialContributions(
678                         enerd->dvdl_lin, lambda, *inputrec.fepvals);
679                 std::tie(foreignLambdaDeltaH, foreignLambdaDhDl) = enerd->foreignLambdaTerms.getTerms(cr);
680             }
681
682             auto& forceWithVirial = (mtsLevel == 0) ? forceWithVirialMtsLevel0 : forceWithVirialMtsLevel1;
683             enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(
684                     inputrec.pbcType,
685                     gmx::arrayRefFromArray(mdatoms->massT, mdatoms->nr),
686                     foreignLambdaDeltaH,
687                     foreignLambdaDhDl,
688                     box,
689                     forceWithVirial,
690                     t,
691                     step,
692                     wcycle,
693                     fplog);
694         }
695     }
696
697     rvec* f = as_rvec_array(forceWithVirialMtsLevel0->force_.data());
698
699     /* Add the forces from enforced rotation potentials (if any) */
700     if (inputrec.bRot)
701     {
702         wallcycle_start(wcycle, ewcROTadd);
703         enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
704         wallcycle_stop(wcycle, ewcROTadd);
705     }
706
707     if (ed)
708     {
709         /* Note that since init_edsam() is called after the initialization
710          * of forcerec, edsam doesn't request the noVirSum force buffer.
711          * Thus if no other algorithm (e.g. PME) requires it, the forces
712          * here will contribute to the virial.
713          */
714         do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
715     }
716
717     /* Add forces from interactive molecular dynamics (IMD), if any */
718     if (inputrec.bIMD && stepWork.computeForces)
719     {
720         imdSession->applyForces(f);
721     }
722 }
723
724 /*! \brief Launch the prepare_step and spread stages of PME GPU.
725  *
726  * \param[in]  pmedata              The PME structure
727  * \param[in]  box                  The box matrix
728  * \param[in]  stepWork             Step schedule flags
729  * \param[in]  xReadyOnDevice       Event synchronizer indicating that the coordinates are ready in the device memory.
730  * \param[in]  lambdaQ              The Coulomb lambda of the current state.
731  * \param[in]  wcycle               The wallcycle structure
732  */
733 static inline void launchPmeGpuSpread(gmx_pme_t*            pmedata,
734                                       const matrix          box,
735                                       const StepWorkload&   stepWork,
736                                       GpuEventSynchronizer* xReadyOnDevice,
737                                       const real            lambdaQ,
738                                       gmx_wallcycle_t       wcycle)
739 {
740     pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
741     pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle, lambdaQ);
742 }
743
744 /*! \brief Launch the FFT and gather stages of PME GPU
745  *
746  * This function only implements setting the output forces (no accumulation).
747  *
748  * \param[in]  pmedata        The PME structure
749  * \param[in]  lambdaQ        The Coulomb lambda of the current system state.
750  * \param[in]  wcycle         The wallcycle structure
751  * \param[in]  stepWork       Step schedule flags
752  */
753 static void launchPmeGpuFftAndGather(gmx_pme_t*               pmedata,
754                                      const real               lambdaQ,
755                                      gmx_wallcycle_t          wcycle,
756                                      const gmx::StepWorkload& stepWork)
757 {
758     pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
759     pme_gpu_launch_gather(pmedata, wcycle, lambdaQ);
760 }
761
762 /*! \brief
763  *  Polling wait for either of the PME or nonbonded GPU tasks.
764  *
765  * Instead of a static order in waiting for GPU tasks, this function
766  * polls checking which of the two tasks completes first, and does the
767  * associated force buffer reduction overlapped with the other task.
768  * By doing that, unlike static scheduling order, it can always overlap
769  * one of the reductions, regardless of the GPU task completion order.
770  *
771  * \param[in]     nbv              Nonbonded verlet structure
772  * \param[in,out] pmedata          PME module data
773  * \param[in,out] forceOutputsNonbonded  Force outputs for the non-bonded forces and shift forces
774  * \param[in,out] forceOutputsPme  Force outputs for the PME forces and virial
775  * \param[in,out] enerd            Energy data structure results are reduced into
776  * \param[in]     lambdaQ          The Coulomb lambda of the current system state.
777  * \param[in]     stepWork         Step schedule flags
778  * \param[in]     wcycle           The wallcycle structure
779  */
780 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
781                                         gmx_pme_t*          pmedata,
782                                         gmx::ForceOutputs*  forceOutputsNonbonded,
783                                         gmx::ForceOutputs*  forceOutputsPme,
784                                         gmx_enerdata_t*     enerd,
785                                         const real          lambdaQ,
786                                         const StepWorkload& stepWork,
787                                         gmx_wallcycle_t     wcycle)
788 {
789     bool isPmeGpuDone = false;
790     bool isNbGpuDone  = false;
791
792     gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
793
794     while (!isPmeGpuDone || !isNbGpuDone)
795     {
796         if (!isPmeGpuDone)
797         {
798             GpuTaskCompletion completionType =
799                     (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
800             isPmeGpuDone = pme_gpu_try_finish_task(
801                     pmedata, stepWork, wcycle, &forceOutputsPme->forceWithVirial(), enerd, lambdaQ, completionType);
802         }
803
804         if (!isNbGpuDone)
805         {
806             auto&             forceBuffersNonbonded = forceOutputsNonbonded->forceWithShiftForces();
807             GpuTaskCompletion completionType =
808                     (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
809             isNbGpuDone = Nbnxm::gpu_try_finish_task(
810                     nbv->gpu_nbv,
811                     stepWork,
812                     AtomLocality::Local,
813                     enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
814                     enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
815                     forceBuffersNonbonded.shiftForces(),
816                     completionType,
817                     wcycle);
818
819             if (isNbGpuDone)
820             {
821                 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceBuffersNonbonded.force());
822             }
823         }
824     }
825 }
826
827 /*! \brief Set up the different force buffers; also does clearing.
828  *
829  * \param[in] forceHelperBuffers  Helper force buffers
830  * \param[in] force     force array
831  * \param[in] stepWork  Step schedule flags
832  * \param[out] wcycle   wallcycle recording structure
833  *
834  * \returns             Cleared force output structure
835  */
836 static ForceOutputs setupForceOutputs(ForceHelperBuffers*                 forceHelperBuffers,
837                                       gmx::ArrayRefWithPadding<gmx::RVec> force,
838                                       const StepWorkload&                 stepWork,
839                                       gmx_wallcycle_t                     wcycle)
840 {
841     wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
842
843     /* NOTE: We assume fr->shiftForces is all zeros here */
844     gmx::ForceWithShiftForces forceWithShiftForces(
845             force, stepWork.computeVirial, forceHelperBuffers->shiftForces());
846
847     if (stepWork.computeForces)
848     {
849         /* Clear the short- and long-range forces */
850         clearRVecs(forceWithShiftForces.force(), true);
851
852         /* Clear the shift forces */
853         clearRVecs(forceWithShiftForces.shiftForces(), false);
854     }
855
856     /* If we need to compute the virial, we might need a separate
857      * force buffer for algorithms for which the virial is calculated
858      * directly, such as PME. Otherwise, forceWithVirial uses the
859      * the same force (f in legacy calls) buffer as other algorithms.
860      */
861     const bool useSeparateForceWithVirialBuffer =
862             (stepWork.computeForces
863              && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
864     /* forceWithVirial uses the local atom range only */
865     gmx::ForceWithVirial forceWithVirial(
866             useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
867                                              : force.unpaddedArrayRef(),
868             stepWork.computeVirial);
869
870     if (useSeparateForceWithVirialBuffer)
871     {
872         /* TODO: update comment
873          * We only compute forces on local atoms. Note that vsites can
874          * spread to non-local atoms, but that part of the buffer is
875          * cleared separately in the vsite spreading code.
876          */
877         clearRVecs(forceWithVirial.force_, true);
878     }
879
880     wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
881
882     return ForceOutputs(
883             forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(), forceWithVirial);
884 }
885
886
887 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
888  */
889 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec&         inputrec,
890                                                           const t_forcerec&         fr,
891                                                           const pull_t*             pull_work,
892                                                           const gmx_edsam*          ed,
893                                                           const t_mdatoms&          mdatoms,
894                                                           const SimulationWorkload& simulationWork,
895                                                           const StepWorkload&       stepWork)
896 {
897     DomainLifetimeWorkload domainWork;
898     // Note that haveSpecialForces is constant over the whole run
899     domainWork.haveSpecialForces =
900             haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
901     domainWork.haveCpuListedForceWork = false;
902     domainWork.haveCpuBondedWork      = false;
903     for (const auto& listedForces : fr.listedForces)
904     {
905         if (listedForces.haveCpuListedForces(*fr.fcdata))
906         {
907             domainWork.haveCpuListedForceWork = true;
908         }
909         if (listedForces.haveCpuBondeds())
910         {
911             domainWork.haveCpuBondedWork = true;
912         }
913     }
914     domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
915     // Note that haveFreeEnergyWork is constant over the whole run
916     domainWork.haveFreeEnergyWork =
917             (fr.efep != FreeEnergyPerturbationType::No && mdatoms.nPerturbed != 0);
918     // We assume we have local force work if there are CPU
919     // force tasks including PME or nonbondeds.
920     domainWork.haveCpuLocalForceWork =
921             domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
922             || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
923             || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
924
925     return domainWork;
926 }
927
928 /*! \brief Set up force flag stuct from the force bitmask.
929  *
930  * \param[in]      legacyFlags          Force bitmask flags used to construct the new flags
931  * \param[in]      mtsLevels            The multiple time-stepping levels, either empty or 2 levels
932  * \param[in]      step                 The current MD step
933  * \param[in]      simulationWork       Simulation workload description.
934  * \param[in]      rankHasPmeDuty       If this rank computes PME.
935  *
936  * \returns New Stepworkload description.
937  */
938 static StepWorkload setupStepWorkload(const int                     legacyFlags,
939                                       ArrayRef<const gmx::MtsLevel> mtsLevels,
940                                       const int64_t                 step,
941                                       const SimulationWorkload&     simulationWork,
942                                       const bool                    rankHasPmeDuty)
943 {
944     GMX_ASSERT(mtsLevels.empty() || mtsLevels.size() == 2, "Expect 0 or 2 MTS levels");
945     const bool computeSlowForces = (mtsLevels.empty() || step % mtsLevels[1].stepFactor == 0);
946
947     StepWorkload flags;
948     flags.stateChanged        = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
949     flags.haveDynamicBox      = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
950     flags.doNeighborSearch    = ((legacyFlags & GMX_FORCE_NS) != 0);
951     flags.computeSlowForces   = computeSlowForces;
952     flags.computeVirial       = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
953     flags.computeEnergy       = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
954     flags.computeForces       = ((legacyFlags & GMX_FORCE_FORCES) != 0);
955     flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
956     flags.computeNonbondedForces =
957             ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded
958             && !(simulationWork.computeNonbondedAtMtsLevel1 && !computeSlowForces);
959     flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
960
961     if (simulationWork.useGpuBufferOps)
962     {
963         GMX_ASSERT(simulationWork.useGpuNonbonded,
964                    "Can only offload buffer ops if nonbonded computation is also offloaded");
965     }
966     flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
967     // on virial steps the CPU reduction path is taken
968     flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
969     flags.useGpuPmeFReduction = flags.computeSlowForces && flags.useGpuFBufferOps && simulationWork.useGpuPme
970                                 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication);
971     flags.useGpuXHalo = simulationWork.useGpuHaloExchange;
972     flags.useGpuFHalo = simulationWork.useGpuHaloExchange && flags.useGpuFBufferOps;
973
974     return flags;
975 }
976
977
978 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
979  *
980  * TODO: eliminate \p useGpuPmeOnThisRank when this is
981  * incorporated in DomainLifetimeWorkload.
982  */
983 static void launchGpuEndOfStepTasks(nonbonded_verlet_t*               nbv,
984                                     gmx::GpuBonded*                   gpuBonded,
985                                     gmx_pme_t*                        pmedata,
986                                     gmx_enerdata_t*                   enerd,
987                                     const gmx::MdrunScheduleWorkload& runScheduleWork,
988                                     bool                              useGpuPmeOnThisRank,
989                                     int64_t                           step,
990                                     gmx_wallcycle_t                   wcycle)
991 {
992     if (runScheduleWork.simulationWork.useGpuNonbonded && runScheduleWork.stepWork.computeNonbondedForces)
993     {
994         /* Launch pruning before buffer clearing because the API overhead of the
995          * clear kernel launches can leave the GPU idle while it could be running
996          * the prune kernel.
997          */
998         if (nbv->isDynamicPruningStepGpu(step))
999         {
1000             nbv->dispatchPruneKernelGpu(step);
1001         }
1002
1003         /* now clear the GPU outputs while we finish the step on the CPU */
1004         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1005         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1006         Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
1007         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1008         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1009     }
1010
1011     if (useGpuPmeOnThisRank)
1012     {
1013         pme_gpu_reinit_computation(pmedata, wcycle);
1014     }
1015
1016     if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
1017     {
1018         // in principle this should be included in the DD balancing region,
1019         // but generally it is infrequent so we'll omit it for the sake of
1020         // simpler code
1021         gpuBonded->waitAccumulateEnergyTerms(enerd);
1022
1023         gpuBonded->clearEnergies();
1024     }
1025 }
1026
1027 //! \brief Data structure to hold dipole-related data and staging arrays
1028 struct DipoleData
1029 {
1030     //! Dipole staging for fast summing over MPI
1031     gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
1032     //! Dipole staging for states A and B (index 0 and 1 resp.)
1033     gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
1034 };
1035
1036
1037 static void reduceAndUpdateMuTot(DipoleData*                   dipoleData,
1038                                  const t_commrec*              cr,
1039                                  const bool                    haveFreeEnergy,
1040                                  gmx::ArrayRef<const real>     lambda,
1041                                  rvec                          muTotal,
1042                                  const DDBalanceRegionHandler& ddBalanceRegionHandler)
1043 {
1044     if (PAR(cr))
1045     {
1046         gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
1047         ddBalanceRegionHandler.reopenRegionCpu();
1048     }
1049     for (int i = 0; i < 2; i++)
1050     {
1051         for (int j = 0; j < DIM; j++)
1052         {
1053             dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
1054         }
1055     }
1056
1057     if (!haveFreeEnergy)
1058     {
1059         copy_rvec(dipoleData->muStateAB[0], muTotal);
1060     }
1061     else
1062     {
1063         for (int j = 0; j < DIM; j++)
1064         {
1065             muTotal[j] = (1.0 - lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)])
1066                                  * dipoleData->muStateAB[0][j]
1067                          + lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]
1068                                    * dipoleData->muStateAB[1][j];
1069         }
1070     }
1071 }
1072
1073 /*! \brief Combines MTS level0 and level1 force buffes into a full and MTS-combined force buffer.
1074  *
1075  * \param[in]     numAtoms        The number of atoms to combine forces for
1076  * \param[in,out] forceMtsLevel0  Input: F_level0, output: F_level0 + F_level1
1077  * \param[in,out] forceMts        Input: F_level1, output: F_level0 + mtsFactor * F_level1
1078  * \param[in]     mtsFactor       The factor between the level0 and level1 time step
1079  */
1080 static void combineMtsForces(const int      numAtoms,
1081                              ArrayRef<RVec> forceMtsLevel0,
1082                              ArrayRef<RVec> forceMts,
1083                              const real     mtsFactor)
1084 {
1085     const int gmx_unused numThreads = gmx_omp_nthreads_get(emntDefault);
1086 #pragma omp parallel for num_threads(numThreads) schedule(static)
1087     for (int i = 0; i < numAtoms; i++)
1088     {
1089         const RVec forceMtsLevel0Tmp = forceMtsLevel0[i];
1090         forceMtsLevel0[i] += forceMts[i];
1091         forceMts[i] = forceMtsLevel0Tmp + mtsFactor * forceMts[i];
1092     }
1093 }
1094
1095 /*! \brief Setup for the local and non-local GPU force reductions:
1096  * reinitialization plus the registration of forces and dependencies.
1097  *
1098  * \param [in] runScheduleWork               Schedule workload flag structure
1099  * \param [in] cr                            Communication record object
1100  * \param [in] fr                            Force record object
1101  */
1102 static void setupGpuForceReductions(gmx::MdrunScheduleWorkload* runScheduleWork,
1103                                     const t_commrec*            cr,
1104                                     t_forcerec*                 fr)
1105 {
1106
1107     nonbonded_verlet_t*          nbv      = fr->nbv.get();
1108     gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1109
1110     // (re-)initialize local GPU force reduction
1111     const bool accumulate =
1112             runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr);
1113     const int atomStart = 0;
1114     fr->gpuForceReduction[gmx::AtomLocality::Local]->reinit(stateGpu->getForces(),
1115                                                             nbv->getNumAtoms(AtomLocality::Local),
1116                                                             nbv->getGridIndices(),
1117                                                             atomStart,
1118                                                             accumulate,
1119                                                             stateGpu->fReducedOnDevice());
1120
1121     // register forces and add dependencies
1122     fr->gpuForceReduction[gmx::AtomLocality::Local]->registerNbnxmForce(nbv->getGpuForces());
1123
1124     if (runScheduleWork->simulationWork.useGpuPme
1125         && (thisRankHasDuty(cr, DUTY_PME) || runScheduleWork->simulationWork.useGpuPmePpCommunication))
1126     {
1127         DeviceBuffer<gmx::RVec> forcePtr =
1128                 thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1129                                               :                    // PME force buffer on same GPU
1130                         fr->pmePpCommGpu->getGpuForceStagingPtr(); // buffer received from other GPU
1131         fr->gpuForceReduction[gmx::AtomLocality::Local]->registerRvecForce(forcePtr);
1132
1133         GpuEventSynchronizer* const pmeSynchronizer =
1134                 (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1135                                                : // PME force buffer on same GPU
1136                          fr->pmePpCommGpu->getForcesReadySynchronizer()); // buffer received from other GPU
1137         fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(pmeSynchronizer);
1138     }
1139
1140     if ((runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr))
1141         && !runScheduleWork->simulationWork.useGpuHaloExchange)
1142     {
1143         auto forcesReadyLocality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
1144         const bool useGpuForceBufferOps = true;
1145         fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1146                 stateGpu->getForcesReadyOnDeviceEvent(forcesReadyLocality, useGpuForceBufferOps));
1147     }
1148
1149     if (runScheduleWork->simulationWork.useGpuHaloExchange)
1150     {
1151         fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1152                 cr->dd->gpuHaloExchange[0][0]->getForcesReadyOnDeviceEvent());
1153     }
1154
1155     if (havePPDomainDecomposition(cr))
1156     {
1157         // (re-)initialize non-local GPU force reduction
1158         const bool accumulate = runScheduleWork->domainWork.haveCpuBondedWork
1159                                 || runScheduleWork->domainWork.haveFreeEnergyWork;
1160         const int atomStart = dd_numHomeAtoms(*cr->dd);
1161         fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->reinit(stateGpu->getForces(),
1162                                                                    nbv->getNumAtoms(AtomLocality::NonLocal),
1163                                                                    nbv->getGridIndices(),
1164                                                                    atomStart,
1165                                                                    accumulate);
1166
1167         // register forces and add dependencies
1168         fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->registerNbnxmForce(nbv->getGpuForces());
1169         if (runScheduleWork->domainWork.haveCpuBondedWork || runScheduleWork->domainWork.haveFreeEnergyWork)
1170         {
1171             fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->addDependency(
1172                     stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::NonLocal, true));
1173         }
1174     }
1175 }
1176
1177
1178 void do_force(FILE*                               fplog,
1179               const t_commrec*                    cr,
1180               const gmx_multisim_t*               ms,
1181               const t_inputrec&                   inputrec,
1182               gmx::Awh*                           awh,
1183               gmx_enfrot*                         enforcedRotation,
1184               gmx::ImdSession*                    imdSession,
1185               pull_t*                             pull_work,
1186               int64_t                             step,
1187               t_nrnb*                             nrnb,
1188               gmx_wallcycle_t                     wcycle,
1189               const gmx_localtop_t*               top,
1190               const matrix                        box,
1191               gmx::ArrayRefWithPadding<gmx::RVec> x,
1192               const history_t*                    hist,
1193               gmx::ForceBuffersView*              forceView,
1194               tensor                              vir_force,
1195               const t_mdatoms*                    mdatoms,
1196               gmx_enerdata_t*                     enerd,
1197               gmx::ArrayRef<const real>           lambda,
1198               t_forcerec*                         fr,
1199               gmx::MdrunScheduleWorkload*         runScheduleWork,
1200               gmx::VirtualSitesHandler*           vsite,
1201               rvec                                muTotal,
1202               double                              t,
1203               gmx_edsam*                          ed,
1204               int                                 legacyFlags,
1205               const DDBalanceRegionHandler&       ddBalanceRegionHandler)
1206 {
1207     auto force = forceView->forceWithPadding();
1208     GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1209                "The size of the force buffer should be at least the number of atoms to compute "
1210                "forces for");
1211
1212     nonbonded_verlet_t*  nbv = fr->nbv.get();
1213     interaction_const_t* ic  = fr->ic.get();
1214
1215     gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1216
1217     const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1218
1219     runScheduleWork->stepWork = setupStepWorkload(
1220             legacyFlags, inputrec.mtsLevels, step, simulationWork, thisRankHasDuty(cr, DUTY_PME));
1221     const StepWorkload& stepWork = runScheduleWork->stepWork;
1222
1223     const bool useGpuPmeOnThisRank =
1224             simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces;
1225
1226     /* At a search step we need to start the first balancing region
1227      * somewhere early inside the step after communication during domain
1228      * decomposition (and not during the previous step as usual).
1229      */
1230     if (stepWork.doNeighborSearch)
1231     {
1232         ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1233     }
1234
1235     clear_mat(vir_force);
1236
1237     if (fr->pbcType != PbcType::No)
1238     {
1239         /* Compute shift vectors every step,
1240          * because of pressure coupling or box deformation!
1241          */
1242         if (stepWork.haveDynamicBox && stepWork.stateChanged)
1243         {
1244             calc_shifts(box, fr->shift_vec);
1245         }
1246
1247         const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1248         const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1249         if (calcCGCM)
1250         {
1251             put_atoms_in_box_omp(fr->pbcType,
1252                                  box,
1253                                  x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1254                                  gmx_omp_nthreads_get(emntDefault));
1255             inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1256         }
1257     }
1258
1259     nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1260
1261     const bool pmeSendCoordinatesFromGpu =
1262             GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1263     const bool reinitGpuPmePpComms =
1264             GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1265
1266     const auto localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1267                                              ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1268                                                        AtomLocality::Local, simulationWork, stepWork)
1269                                              : nullptr;
1270
1271     // Copy coordinate from the GPU if update is on the GPU and there
1272     // are forces to be computed on the CPU, or for the computation of
1273     // virial, or if host-side data will be transferred from this task
1274     // to a remote task for halo exchange or PME-PP communication. At
1275     // search steps the current coordinates are already on the host,
1276     // hence copy is not needed.
1277     const bool haveHostPmePpComms =
1278             !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1279
1280     GMX_ASSERT(simulationWork.useGpuHaloExchange
1281                        == ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange[0].empty())),
1282                "The GPU halo exchange is active, but it has not been constructed.");
1283     const bool haveHostHaloExchangeComms =
1284             havePPDomainDecomposition(cr) && !simulationWork.useGpuHaloExchange;
1285
1286     bool gmx_used_in_debug haveCopiedXFromGpu = false;
1287     if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1288         && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1289             || haveHostPmePpComms || haveHostHaloExchangeComms))
1290     {
1291         stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1292         haveCopiedXFromGpu = true;
1293     }
1294
1295     // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1296     // Otherwise the send will occur after H2D coordinate transfer.
1297     if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu && stepWork.computeSlowForces)
1298     {
1299         /* Send particle coordinates to the pme nodes */
1300         if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1301         {
1302             GMX_ASSERT(haveCopiedXFromGpu,
1303                        "a wait should only be triggered if copy has been scheduled");
1304             stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1305         }
1306
1307         gmx_pme_send_coordinates(fr,
1308                                  cr,
1309                                  box,
1310                                  as_rvec_array(x.unpaddedArrayRef().data()),
1311                                  lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1312                                  lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1313                                  (stepWork.computeVirial || stepWork.computeEnergy),
1314                                  step,
1315                                  simulationWork.useGpuPmePpCommunication,
1316                                  reinitGpuPmePpComms,
1317                                  pmeSendCoordinatesFromGpu,
1318                                  localXReadyOnDevice,
1319                                  wcycle);
1320     }
1321
1322     // Coordinates on the device are needed if PME or BufferOps are offloaded.
1323     // The local coordinates can be copied right away.
1324     // NOTE: Consider moving this copy to right after they are updated and constrained,
1325     //       if the later is not offloaded.
1326     if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1327     {
1328         if (stepWork.doNeighborSearch)
1329         {
1330             // TODO refactor this to do_md, after partitioning.
1331             stateGpu->reinit(mdatoms->homenr,
1332                              cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1333             if (useGpuPmeOnThisRank)
1334             {
1335                 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1336                 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1337             }
1338         }
1339         // We need to copy coordinates when:
1340         // 1. Update is not offloaded
1341         // 2. The buffers were reinitialized on search step
1342         if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1343         {
1344             GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1345             stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1346         }
1347     }
1348
1349     // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1350     // Otherwise the send will occur before the H2D coordinate transfer.
1351     if (!thisRankHasDuty(cr, DUTY_PME) && pmeSendCoordinatesFromGpu)
1352     {
1353         /* Send particle coordinates to the pme nodes */
1354         gmx_pme_send_coordinates(fr,
1355                                  cr,
1356                                  box,
1357                                  as_rvec_array(x.unpaddedArrayRef().data()),
1358                                  lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1359                                  lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1360                                  (stepWork.computeVirial || stepWork.computeEnergy),
1361                                  step,
1362                                  simulationWork.useGpuPmePpCommunication,
1363                                  reinitGpuPmePpComms,
1364                                  pmeSendCoordinatesFromGpu,
1365                                  localXReadyOnDevice,
1366                                  wcycle);
1367     }
1368
1369     if (useGpuPmeOnThisRank)
1370     {
1371         launchPmeGpuSpread(fr->pmedata,
1372                            box,
1373                            stepWork,
1374                            localXReadyOnDevice,
1375                            lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1376                            wcycle);
1377     }
1378
1379     const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1380
1381     /* do gridding for pair search */
1382     if (stepWork.doNeighborSearch)
1383     {
1384         if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1385         {
1386             fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1387         }
1388
1389         wallcycle_start(wcycle, ewcNS);
1390         if (!DOMAINDECOMP(cr))
1391         {
1392             const rvec vzero       = { 0.0_real, 0.0_real, 0.0_real };
1393             const rvec boxDiagonal = { box[XX][XX], box[YY][YY], box[ZZ][ZZ] };
1394             wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1395             nbnxn_put_on_grid(nbv,
1396                               box,
1397                               0,
1398                               vzero,
1399                               boxDiagonal,
1400                               nullptr,
1401                               { 0, mdatoms->homenr },
1402                               -1,
1403                               fr->cginfo,
1404                               x.unpaddedArrayRef(),
1405                               0,
1406                               nullptr);
1407             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1408         }
1409         else
1410         {
1411             wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1412             nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1413             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1414         }
1415
1416         nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1417                                gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1418                                fr->cginfo);
1419
1420         wallcycle_stop(wcycle, ewcNS);
1421
1422         /* initialize the GPU nbnxm atom data and bonded data structures */
1423         if (simulationWork.useGpuNonbonded)
1424         {
1425             // Note: cycle counting only nononbondeds, gpuBonded counts internally
1426             wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1427             wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1428             Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1429             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1430             wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1431
1432             if (fr->gpuBonded)
1433             {
1434                 /* Now we put all atoms on the grid, we can assign bonded
1435                  * interactions to the GPU, where the grid order is
1436                  * needed. Also the xq, f and fshift device buffers have
1437                  * been reallocated if needed, so the bonded code can
1438                  * learn about them. */
1439                 // TODO the xq, f, and fshift buffers are now shared
1440                 // resources, so they should be maintained by a
1441                 // higher-level object than the nb module.
1442                 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1443                                                                       top->idef,
1444                                                                       Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1445                                                                       Nbnxm::gpu_get_f(nbv->gpu_nbv),
1446                                                                       Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1447             }
1448         }
1449
1450         // Need to run after the GPU-offload bonded interaction lists
1451         // are set up to be able to determine whether there is bonded work.
1452         runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1453                 inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1454
1455         wallcycle_start_nocount(wcycle, ewcNS);
1456         wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1457         /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1458         nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1459
1460         nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1461
1462         wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1463         wallcycle_stop(wcycle, ewcNS);
1464
1465         if (stepWork.useGpuXBufferOps)
1466         {
1467             nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1468         }
1469
1470         if (simulationWork.useGpuBufferOps)
1471         {
1472             setupGpuForceReductions(runScheduleWork, cr, fr);
1473         }
1474     }
1475     else if (!EI_TPI(inputrec.eI) && stepWork.computeNonbondedForces)
1476     {
1477         if (stepWork.useGpuXBufferOps)
1478         {
1479             GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1480             nbv->convertCoordinatesGpu(AtomLocality::Local, stateGpu->getCoordinates(), localXReadyOnDevice);
1481         }
1482         else
1483         {
1484             if (simulationWork.useGpuUpdate)
1485             {
1486                 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1487                 GMX_ASSERT(haveCopiedXFromGpu,
1488                            "a wait should only be triggered if copy has been scheduled");
1489                 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1490             }
1491             nbv->convertCoordinates(AtomLocality::Local, x.unpaddedArrayRef());
1492         }
1493     }
1494
1495     if (simulationWork.useGpuNonbonded && (stepWork.computeNonbondedForces || domainWork.haveGpuBondedWork))
1496     {
1497         ddBalanceRegionHandler.openBeforeForceComputationGpu();
1498
1499         wallcycle_start(wcycle, ewcLAUNCH_GPU);
1500         wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1501         Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1502         if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1503         {
1504             Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1505         }
1506         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1507         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1508         // with X buffer ops offloaded to the GPU on all but the search steps
1509
1510         // bonded work not split into separate local and non-local, so with DD
1511         // we can only launch the kernel after non-local coordinates have been received.
1512         if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1513         {
1514             fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1515         }
1516
1517         /* launch local nonbonded work on GPU */
1518         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1519         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1520         do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1521         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1522         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1523     }
1524
1525     if (useGpuPmeOnThisRank)
1526     {
1527         // In PME GPU and mixed mode we launch FFT / gather after the
1528         // X copy/transform to allow overlap as well as after the GPU NB
1529         // launch to avoid FFT launch overhead hijacking the CPU and delaying
1530         // the nonbonded kernel.
1531         launchPmeGpuFftAndGather(fr->pmedata,
1532                                  lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1533                                  wcycle,
1534                                  stepWork);
1535     }
1536
1537     /* Communicate coordinates and sum dipole if necessary +
1538        do non-local pair search */
1539     if (havePPDomainDecomposition(cr))
1540     {
1541         if (stepWork.doNeighborSearch)
1542         {
1543             // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1544             wallcycle_start_nocount(wcycle, ewcNS);
1545             wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1546             /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1547             nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1548
1549             nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1550             wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1551             wallcycle_stop(wcycle, ewcNS);
1552             // TODO refactor this GPU halo exchange re-initialisation
1553             // to location in do_md where GPU halo exchange is
1554             // constructed at partitioning, after above stateGpu
1555             // re-initialization has similarly been refactored
1556             if (simulationWork.useGpuHaloExchange)
1557             {
1558                 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1559             }
1560         }
1561         else
1562         {
1563             if (stepWork.useGpuXHalo)
1564             {
1565                 // The following must be called after local setCoordinates (which records an event
1566                 // when the coordinate data has been copied to the device).
1567                 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1568
1569                 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1570                 {
1571                     // non-local part of coordinate buffer must be copied back to host for CPU work
1572                     stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1573                 }
1574             }
1575             else
1576             {
1577                 if (simulationWork.useGpuUpdate)
1578                 {
1579                     GMX_ASSERT(haveCopiedXFromGpu,
1580                                "a wait should only be triggered if copy has been scheduled");
1581                     stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1582                 }
1583                 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1584             }
1585
1586             if (stepWork.useGpuXBufferOps)
1587             {
1588                 if (!useGpuPmeOnThisRank && !stepWork.useGpuXHalo)
1589                 {
1590                     stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1591                 }
1592                 nbv->convertCoordinatesGpu(AtomLocality::NonLocal,
1593                                            stateGpu->getCoordinates(),
1594                                            stateGpu->getCoordinatesReadyOnDeviceEvent(
1595                                                    AtomLocality::NonLocal, simulationWork, stepWork));
1596             }
1597             else
1598             {
1599                 nbv->convertCoordinates(AtomLocality::NonLocal, x.unpaddedArrayRef());
1600             }
1601         }
1602
1603         if (simulationWork.useGpuNonbonded)
1604         {
1605
1606             if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1607             {
1608                 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1609                 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1610                 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1611                 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1612                 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1613             }
1614
1615             if (domainWork.haveGpuBondedWork)
1616             {
1617                 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1618             }
1619
1620             /* launch non-local nonbonded tasks on GPU */
1621             wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1622             wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1623             do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1624             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1625             wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1626         }
1627     }
1628
1629     if (simulationWork.useGpuNonbonded && stepWork.computeNonbondedForces)
1630     {
1631         /* launch D2H copy-back F */
1632         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1633         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1634
1635         if (havePPDomainDecomposition(cr))
1636         {
1637             Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1638         }
1639         Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1640         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1641
1642         if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1643         {
1644             fr->gpuBonded->launchEnergyTransfer();
1645         }
1646         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1647     }
1648
1649     gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1650     if (fr->wholeMoleculeTransform)
1651     {
1652         xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1653     }
1654
1655     DipoleData dipoleData;
1656
1657     if (simulationWork.computeMuTot)
1658     {
1659         const int start = 0;
1660
1661         /* Calculate total (local) dipole moment in a temporary common array.
1662          * This makes it possible to sum them over nodes faster.
1663          */
1664         gmx::ArrayRef<const gmx::RVec> xRef =
1665                 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1666         calc_mu(start,
1667                 mdatoms->homenr,
1668                 xRef,
1669                 gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1670                 gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
1671                 mdatoms->nChargePerturbed != 0,
1672                 dipoleData.muStaging[0],
1673                 dipoleData.muStaging[1]);
1674
1675         reduceAndUpdateMuTot(
1676                 &dipoleData, cr, (fr->efep != FreeEnergyPerturbationType::No), lambda, muTotal, ddBalanceRegionHandler);
1677     }
1678
1679     /* Reset energies */
1680     reset_enerdata(enerd);
1681
1682     if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1683     {
1684         wallcycle_start(wcycle, ewcPPDURINGPME);
1685         dd_force_flop_start(cr->dd, nrnb);
1686     }
1687
1688     // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1689     // this wait ensures that the D2H transfer is complete.
1690     if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1691         && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1692     {
1693         GMX_ASSERT(haveCopiedXFromGpu, "a wait should only be triggered if copy has been scheduled");
1694         stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1695     }
1696
1697     if (inputrec.bRot)
1698     {
1699         wallcycle_start(wcycle, ewcROT);
1700         do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, stepWork.doNeighborSearch);
1701         wallcycle_stop(wcycle, ewcROT);
1702     }
1703
1704     /* Start the force cycle counter.
1705      * Note that a different counter is used for dynamic load balancing.
1706      */
1707     wallcycle_start(wcycle, ewcFORCE);
1708
1709     /* Set up and clear force outputs:
1710      * forceOutMtsLevel0:  everything except what is in the other two outputs
1711      * forceOutMtsLevel1:  PME-mesh and listed-forces group 1
1712      * forceOutNonbonded: non-bonded forces
1713      * Without multiple time stepping all point to the same object.
1714      * With multiple time-stepping the use is different for MTS fast (level0 only) and slow steps.
1715      */
1716     ForceOutputs forceOutMtsLevel0 =
1717             setupForceOutputs(&fr->forceHelperBuffers[0], force, stepWork, wcycle);
1718
1719     // Force output for MTS combined forces, only set at level1 MTS steps
1720     std::optional<ForceOutputs> forceOutMts =
1721             (fr->useMts && stepWork.computeSlowForces)
1722                     ? std::optional(setupForceOutputs(&fr->forceHelperBuffers[1],
1723                                                       forceView->forceMtsCombinedWithPadding(),
1724                                                       stepWork,
1725                                                       wcycle))
1726                     : std::nullopt;
1727
1728     ForceOutputs* forceOutMtsLevel1 =
1729             fr->useMts ? (stepWork.computeSlowForces ? &forceOutMts.value() : nullptr) : &forceOutMtsLevel0;
1730
1731     const bool nonbondedAtMtsLevel1 = runScheduleWork->simulationWork.computeNonbondedAtMtsLevel1;
1732
1733     ForceOutputs* forceOutNonbonded = nonbondedAtMtsLevel1 ? forceOutMtsLevel1 : &forceOutMtsLevel0;
1734
1735     if (inputrec.bPull && pull_have_constraint(*pull_work))
1736     {
1737         clear_pull_forces(pull_work);
1738     }
1739
1740     /* We calculate the non-bonded forces, when done on the CPU, here.
1741      * We do this before calling do_force_lowlevel, because in that
1742      * function, the listed forces are calculated before PME, which
1743      * does communication.  With this order, non-bonded and listed
1744      * force calculation imbalance can be balanced out by the domain
1745      * decomposition load balancing.
1746      */
1747
1748     const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1749
1750     if (!useOrEmulateGpuNb)
1751     {
1752         do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1753     }
1754
1755     if (fr->efep != FreeEnergyPerturbationType::No && stepWork.computeNonbondedForces)
1756     {
1757         /* Calculate the local and non-local free energy interactions here.
1758          * Happens here on the CPU both with and without GPU.
1759          */
1760         nbv->dispatchFreeEnergyKernel(InteractionLocality::Local,
1761                                       fr,
1762                                       as_rvec_array(x.unpaddedArrayRef().data()),
1763                                       &forceOutNonbonded->forceWithShiftForces(),
1764                                       *mdatoms,
1765                                       inputrec.fepvals.get(),
1766                                       lambda,
1767                                       enerd,
1768                                       stepWork,
1769                                       nrnb);
1770
1771         if (havePPDomainDecomposition(cr))
1772         {
1773             nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal,
1774                                           fr,
1775                                           as_rvec_array(x.unpaddedArrayRef().data()),
1776                                           &forceOutNonbonded->forceWithShiftForces(),
1777                                           *mdatoms,
1778                                           inputrec.fepvals.get(),
1779                                           lambda,
1780                                           enerd,
1781                                           stepWork,
1782                                           nrnb);
1783         }
1784     }
1785
1786     if (stepWork.computeNonbondedForces && !useOrEmulateGpuNb)
1787     {
1788         if (havePPDomainDecomposition(cr))
1789         {
1790             do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1791         }
1792
1793         if (stepWork.computeForces)
1794         {
1795             /* Add all the non-bonded force to the normal force array.
1796              * This can be split into a local and a non-local part when overlapping
1797              * communication with calculation with domain decomposition.
1798              */
1799             wallcycle_stop(wcycle, ewcFORCE);
1800             nbv->atomdata_add_nbat_f_to_f(AtomLocality::All,
1801                                           forceOutNonbonded->forceWithShiftForces().force());
1802             wallcycle_start_nocount(wcycle, ewcFORCE);
1803         }
1804
1805         /* If there are multiple fshift output buffers we need to reduce them */
1806         if (stepWork.computeVirial)
1807         {
1808             /* This is not in a subcounter because it takes a
1809                negligible and constant-sized amount of time */
1810             nbnxn_atomdata_add_nbat_fshift_to_fshift(
1811                     *nbv->nbat, forceOutNonbonded->forceWithShiftForces().shiftForces());
1812         }
1813     }
1814
1815     // TODO Force flags should include haveFreeEnergyWork for this domain
1816     if (stepWork.useGpuXHalo && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1817     {
1818         wallcycle_stop(wcycle, ewcFORCE);
1819         /* Wait for non-local coordinate data to be copied from device */
1820         stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1821         wallcycle_start_nocount(wcycle, ewcFORCE);
1822     }
1823
1824     // Compute wall interactions, when present.
1825     // Note: should be moved to special forces.
1826     if (inputrec.nwall && stepWork.computeNonbondedForces)
1827     {
1828         /* foreign lambda component for walls */
1829         real dvdl_walls = do_walls(inputrec,
1830                                    *fr,
1831                                    box,
1832                                    *mdatoms,
1833                                    x.unpaddedConstArrayRef(),
1834                                    &forceOutMtsLevel0.forceWithVirial(),
1835                                    lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1836                                    enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
1837                                    nrnb);
1838         enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += dvdl_walls;
1839     }
1840
1841     if (stepWork.computeListedForces)
1842     {
1843         /* Check whether we need to take into account PBC in listed interactions */
1844         bool needMolPbc = false;
1845         for (const auto& listedForces : fr->listedForces)
1846         {
1847             if (listedForces.haveCpuListedForces(*fr->fcdata))
1848             {
1849                 needMolPbc = fr->bMolPBC;
1850             }
1851         }
1852
1853         t_pbc pbc;
1854
1855         if (needMolPbc)
1856         {
1857             /* Since all atoms are in the rectangular or triclinic unit-cell,
1858              * only single box vector shifts (2 in x) are required.
1859              */
1860             set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1861         }
1862
1863         for (int mtsIndex = 0; mtsIndex < (fr->useMts && stepWork.computeSlowForces ? 2 : 1); mtsIndex++)
1864         {
1865             ListedForces& listedForces = fr->listedForces[mtsIndex];
1866             ForceOutputs& forceOut     = (mtsIndex == 0 ? forceOutMtsLevel0 : *forceOutMtsLevel1);
1867             listedForces.calculate(wcycle,
1868                                    box,
1869                                    inputrec.fepvals.get(),
1870                                    cr,
1871                                    ms,
1872                                    x,
1873                                    xWholeMolecules,
1874                                    fr->fcdata.get(),
1875                                    hist,
1876                                    &forceOut,
1877                                    fr,
1878                                    &pbc,
1879                                    enerd,
1880                                    nrnb,
1881                                    lambda,
1882                                    mdatoms,
1883                                    DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr,
1884                                    stepWork);
1885         }
1886     }
1887
1888     if (stepWork.computeSlowForces)
1889     {
1890         calculateLongRangeNonbondeds(fr,
1891                                      inputrec,
1892                                      cr,
1893                                      nrnb,
1894                                      wcycle,
1895                                      mdatoms,
1896                                      x.unpaddedConstArrayRef(),
1897                                      &forceOutMtsLevel1->forceWithVirial(),
1898                                      enerd,
1899                                      box,
1900                                      lambda,
1901                                      dipoleData.muStateAB,
1902                                      stepWork,
1903                                      ddBalanceRegionHandler);
1904     }
1905
1906     wallcycle_stop(wcycle, ewcFORCE);
1907
1908     // VdW dispersion correction, only computed on master rank to avoid double counting
1909     if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1910     {
1911         // Calculate long range corrections to pressure and energy
1912         const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
1913                 box, lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)]);
1914
1915         if (stepWork.computeEnergy)
1916         {
1917             enerd->term[F_DISPCORR] = correction.energy;
1918             enerd->term[F_DVDL_VDW] += correction.dvdl;
1919             enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += correction.dvdl;
1920         }
1921         if (stepWork.computeVirial)
1922         {
1923             correction.correctVirial(vir_force);
1924             enerd->term[F_PDISPCORR] = correction.pressure;
1925         }
1926     }
1927
1928     computeSpecialForces(fplog,
1929                          cr,
1930                          inputrec,
1931                          awh,
1932                          enforcedRotation,
1933                          imdSession,
1934                          pull_work,
1935                          step,
1936                          t,
1937                          wcycle,
1938                          fr->forceProviders,
1939                          box,
1940                          x.unpaddedArrayRef(),
1941                          mdatoms,
1942                          lambda,
1943                          stepWork,
1944                          &forceOutMtsLevel0.forceWithVirial(),
1945                          forceOutMtsLevel1 ? &forceOutMtsLevel1->forceWithVirial() : nullptr,
1946                          enerd,
1947                          ed,
1948                          stepWork.doNeighborSearch);
1949
1950     if (havePPDomainDecomposition(cr) && stepWork.computeForces && stepWork.useGpuFHalo
1951         && domainWork.haveCpuLocalForceWork)
1952     {
1953         stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(), AtomLocality::Local);
1954     }
1955
1956     GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
1957                "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
1958     GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFHalo),
1959                "The schedule below does not allow for nonbonded MTS with GPU halo exchange");
1960     // Will store the amount of cycles spent waiting for the GPU that
1961     // will be later used in the DLB accounting.
1962     float cycles_wait_gpu = 0;
1963     if (useOrEmulateGpuNb && stepWork.computeNonbondedForces)
1964     {
1965         auto& forceWithShiftForces = forceOutNonbonded->forceWithShiftForces();
1966
1967         /* wait for non-local forces (or calculate in emulation mode) */
1968         if (havePPDomainDecomposition(cr))
1969         {
1970             if (simulationWork.useGpuNonbonded)
1971             {
1972                 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1973                         nbv->gpu_nbv,
1974                         stepWork,
1975                         AtomLocality::NonLocal,
1976                         enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
1977                         enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
1978                         forceWithShiftForces.shiftForces(),
1979                         wcycle);
1980             }
1981             else
1982             {
1983                 wallcycle_start_nocount(wcycle, ewcFORCE);
1984                 do_nb_verlet(
1985                         fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes, step, nrnb, wcycle);
1986                 wallcycle_stop(wcycle, ewcFORCE);
1987             }
1988
1989             if (stepWork.useGpuFBufferOps)
1990             {
1991                 // TODO: move this into DomainLifetimeWorkload, including the second part of the
1992                 // condition The bonded and free energy CPU tasks can have non-local force
1993                 // contributions which are a dependency for the GPU force reduction.
1994                 bool haveNonLocalForceContribInCpuBuffer =
1995                         domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1996
1997                 if (haveNonLocalForceContribInCpuBuffer)
1998                 {
1999                     stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2000                                               AtomLocality::NonLocal);
2001                 }
2002
2003
2004                 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->execute();
2005
2006                 if (!stepWork.useGpuFHalo)
2007                 {
2008                     // copy from GPU input for dd_move_f()
2009                     stateGpu->copyForcesFromGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2010                                                 AtomLocality::NonLocal);
2011                 }
2012             }
2013             else
2014             {
2015                 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
2016             }
2017
2018             if (fr->nbv->emulateGpu() && stepWork.computeVirial)
2019             {
2020                 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
2021             }
2022         }
2023     }
2024
2025     /* Combining the forces for multiple time stepping before the halo exchange, when possible,
2026      * avoids an extra halo exchange (when DD is used) and post-processing step.
2027      */
2028     const bool combineMtsForcesBeforeHaloExchange =
2029             (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2030              && (legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0
2031              && !(stepWork.computeVirial || simulationWork.useGpuNonbonded || useGpuPmeOnThisRank));
2032     if (combineMtsForcesBeforeHaloExchange)
2033     {
2034         const int numAtoms = havePPDomainDecomposition(cr) ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr;
2035         combineMtsForces(numAtoms,
2036                          force.unpaddedArrayRef(),
2037                          forceView->forceMtsCombined(),
2038                          inputrec.mtsLevels[1].stepFactor);
2039     }
2040
2041     if (havePPDomainDecomposition(cr))
2042     {
2043         /* We are done with the CPU compute.
2044          * We will now communicate the non-local forces.
2045          * If we use a GPU this will overlap with GPU work, so in that case
2046          * we do not close the DD force balancing region here.
2047          */
2048         ddBalanceRegionHandler.closeAfterForceComputationCpu();
2049
2050         if (stepWork.computeForces)
2051         {
2052
2053             if (stepWork.useGpuFHalo)
2054             {
2055                 // If there exist CPU forces, data from halo exchange should accumulate into these
2056                 bool accumulateForces = domainWork.haveCpuLocalForceWork;
2057                 if (!accumulateForces)
2058                 {
2059                     // Force halo exchange will set a subset of local atoms with remote non-local data
2060                     // First clear local portion of force array, so that untouched atoms are zero
2061                     stateGpu->clearForcesOnGpu(AtomLocality::Local);
2062                 }
2063                 communicateGpuHaloForces(*cr, accumulateForces);
2064             }
2065             else
2066             {
2067                 if (stepWork.useGpuFBufferOps)
2068                 {
2069                     stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
2070                 }
2071
2072                 // Without MTS or with MTS at slow steps with uncombined forces we need to
2073                 // communicate the fast forces
2074                 if (!fr->useMts || !combineMtsForcesBeforeHaloExchange)
2075                 {
2076                     dd_move_f(cr->dd, &forceOutMtsLevel0.forceWithShiftForces(), wcycle);
2077                 }
2078                 // With MTS we need to communicate the slow or combined (in forceOutMtsLevel1) forces
2079                 if (fr->useMts && stepWork.computeSlowForces)
2080                 {
2081                     dd_move_f(cr->dd, &forceOutMtsLevel1->forceWithShiftForces(), wcycle);
2082                 }
2083             }
2084         }
2085     }
2086
2087     // With both nonbonded and PME offloaded a GPU on the same rank, we use
2088     // an alternating wait/reduction scheme.
2089     bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
2090                              && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
2091     if (alternateGpuWait)
2092     {
2093         alternatePmeNbGpuWaitReduce(fr->nbv.get(),
2094                                     fr->pmedata,
2095                                     forceOutNonbonded,
2096                                     forceOutMtsLevel1,
2097                                     enerd,
2098                                     lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
2099                                     stepWork,
2100                                     wcycle);
2101     }
2102
2103     if (!alternateGpuWait && useGpuPmeOnThisRank)
2104     {
2105         pme_gpu_wait_and_reduce(fr->pmedata,
2106                                 stepWork,
2107                                 wcycle,
2108                                 &forceOutMtsLevel1->forceWithVirial(),
2109                                 enerd,
2110                                 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]);
2111     }
2112
2113     /* Wait for local GPU NB outputs on the non-alternating wait path */
2114     if (!alternateGpuWait && stepWork.computeNonbondedForces && simulationWork.useGpuNonbonded)
2115     {
2116         /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
2117          * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
2118          * but even with a step of 0.1 ms the difference is less than 1%
2119          * of the step time.
2120          */
2121         const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
2122         const float waitCycles               = Nbnxm::gpu_wait_finish_task(
2123                 nbv->gpu_nbv,
2124                 stepWork,
2125                 AtomLocality::Local,
2126                 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
2127                 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
2128                 forceOutNonbonded->forceWithShiftForces().shiftForces(),
2129                 wcycle);
2130
2131         if (ddBalanceRegionHandler.useBalancingRegion())
2132         {
2133             DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
2134             if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
2135             {
2136                 /* We measured few cycles, it could be that the kernel
2137                  * and transfer finished earlier and there was no actual
2138                  * wait time, only API call overhead.
2139                  * Then the actual time could be anywhere between 0 and
2140                  * cycles_wait_est. We will use half of cycles_wait_est.
2141                  */
2142                 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
2143             }
2144             ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
2145         }
2146     }
2147
2148     if (fr->nbv->emulateGpu())
2149     {
2150         // NOTE: emulation kernel is not included in the balancing region,
2151         // but emulation mode does not target performance anyway
2152         wallcycle_start_nocount(wcycle, ewcFORCE);
2153         do_nb_verlet(fr,
2154                      ic,
2155                      enerd,
2156                      stepWork,
2157                      InteractionLocality::Local,
2158                      DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
2159                      step,
2160                      nrnb,
2161                      wcycle);
2162         wallcycle_stop(wcycle, ewcFORCE);
2163     }
2164
2165     // If on GPU PME-PP comms path, receive forces from PME before GPU buffer ops
2166     // TODO refactor this and unify with below default-path call to the same function
2167     if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces
2168         && simulationWork.useGpuPmePpCommunication)
2169     {
2170         /* In case of node-splitting, the PP nodes receive the long-range
2171          * forces, virial and energy from the PME nodes here.
2172          */
2173         pme_receive_force_ener(fr,
2174                                cr,
2175                                &forceOutMtsLevel1->forceWithVirial(),
2176                                enerd,
2177                                simulationWork.useGpuPmePpCommunication,
2178                                stepWork.useGpuPmeFReduction,
2179                                wcycle);
2180     }
2181
2182
2183     /* Do the nonbonded GPU (or emulation) force buffer reduction
2184      * on the non-alternating path. */
2185     GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
2186                "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
2187     if (useOrEmulateGpuNb && !alternateGpuWait)
2188     {
2189         if (stepWork.useGpuFBufferOps)
2190         {
2191             ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2192
2193             // Flag to specify whether the CPU force buffer has contributions to
2194             // local atoms. This depends on whether there are CPU-based force tasks
2195             // or when DD is active the halo exchange has resulted in contributions
2196             // from the non-local part.
2197             const bool haveLocalForceContribInCpuBuffer =
2198                     (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
2199
2200             // TODO: move these steps as early as possible:
2201             // - CPU f H2D should be as soon as all CPU-side forces are done
2202             // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
2203             //   before the next CPU task that consumes the forces: vsite spread or update)
2204             // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
2205             //   of the halo exchange. In that case the copy is instead performed above, before the exchange.
2206             //   These should be unified.
2207             if (haveLocalForceContribInCpuBuffer && !stepWork.useGpuFHalo)
2208             {
2209                 // Note: AtomLocality::All is used for the non-DD case because, as in this
2210                 // case copyForcesToGpu() uses a separate stream, it allows overlap of
2211                 // CPU force H2D with GPU force tasks on all streams including those in the
2212                 // local stream which would otherwise be implicit dependencies for the
2213                 // transfer and would not overlap.
2214                 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
2215
2216                 stateGpu->copyForcesToGpu(forceWithShift, locality);
2217             }
2218
2219             if (stepWork.computeNonbondedForces)
2220             {
2221                 fr->gpuForceReduction[gmx::AtomLocality::Local]->execute();
2222             }
2223
2224             // Copy forces to host if they are needed for update or if virtual sites are enabled.
2225             // If there are vsites, we need to copy forces every step to spread vsite forces on host.
2226             // TODO: When the output flags will be included in step workload, this copy can be combined with the
2227             //       copy call done in sim_utils(...) for the output.
2228             // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
2229             //       they should not be copied in do_md(...) for the output.
2230             if (!simulationWork.useGpuUpdate
2231                 || (simulationWork.useGpuUpdate && DOMAINDECOMP(cr) && haveHostPmePpComms) || vsite)
2232             {
2233                 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
2234                 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
2235             }
2236         }
2237         else if (stepWork.computeNonbondedForces)
2238         {
2239             ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2240             nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
2241         }
2242     }
2243
2244     launchGpuEndOfStepTasks(
2245             nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork, useGpuPmeOnThisRank, step, wcycle);
2246
2247     if (DOMAINDECOMP(cr))
2248     {
2249         dd_force_flop_stop(cr->dd, nrnb);
2250     }
2251
2252     const bool haveCombinedMtsForces = (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2253                                         && combineMtsForcesBeforeHaloExchange);
2254     if (stepWork.computeForces)
2255     {
2256         postProcessForceWithShiftForces(
2257                 nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutMtsLevel0, vir_force, *mdatoms, *fr, vsite, stepWork);
2258
2259         if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2260         {
2261             postProcessForceWithShiftForces(
2262                     nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, *mdatoms, *fr, vsite, stepWork);
2263         }
2264     }
2265
2266     // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
2267     if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
2268         && stepWork.computeSlowForces)
2269     {
2270         /* In case of node-splitting, the PP nodes receive the long-range
2271          * forces, virial and energy from the PME nodes here.
2272          */
2273         pme_receive_force_ener(fr,
2274                                cr,
2275                                &forceOutMtsLevel1->forceWithVirial(),
2276                                enerd,
2277                                simulationWork.useGpuPmePpCommunication,
2278                                false,
2279                                wcycle);
2280     }
2281
2282     if (stepWork.computeForces)
2283     {
2284         /* If we don't use MTS or if we already combined the MTS forces before, we only
2285          * need to post-process one ForceOutputs object here, called forceOutCombined,
2286          * otherwise we have to post-process two outputs and then combine them.
2287          */
2288         ForceOutputs& forceOutCombined = (haveCombinedMtsForces ? forceOutMts.value() : forceOutMtsLevel0);
2289         postProcessForces(
2290                 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutCombined, vir_force, mdatoms, fr, vsite, stepWork);
2291
2292         if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2293         {
2294             postProcessForces(
2295                     cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, mdatoms, fr, vsite, stepWork);
2296
2297             combineMtsForces(mdatoms->homenr,
2298                              force.unpaddedArrayRef(),
2299                              forceView->forceMtsCombined(),
2300                              inputrec.mtsLevels[1].stepFactor);
2301         }
2302     }
2303
2304     if (stepWork.computeEnergy)
2305     {
2306         /* Compute the final potential energy terms */
2307         accumulatePotentialEnergies(enerd, lambda, inputrec.fepvals.get());
2308
2309         if (!EI_TPI(inputrec.eI))
2310         {
2311             checkPotentialEnergyValidity(step, *enerd, inputrec);
2312         }
2313     }
2314
2315     /* In case we don't have constraints and are using GPUs, the next balancing
2316      * region starts here.
2317      * Some "special" work at the end of do_force_cuts?, such as vsite spread,
2318      * virial calculation and COM pulling, is not thus not included in
2319      * the balance timing, which is ok as most tasks do communication.
2320      */
2321     ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);
2322 }