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