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