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