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