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