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