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