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