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