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