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