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