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