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