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