3a827652d996af9d2b8d45fd207a8ecbbebba5c2
[alexxy/gromacs.git] / src / gromacs / mdlib / sim_util.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #include <cmath>
42 #include <cstdint>
43 #include <cstdio>
44 #include <cstring>
45
46 #include <array>
47
48 #include "gromacs/awh/awh.h"
49 #include "gromacs/domdec/dlbtiming.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/domdec/partition.h"
53 #include "gromacs/essentialdynamics/edsam.h"
54 #include "gromacs/ewald/pme.h"
55 #include "gromacs/gmxlib/chargegroup.h"
56 #include "gromacs/gmxlib/network.h"
57 #include "gromacs/gmxlib/nrnb.h"
58 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
59 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
60 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
61 #include "gromacs/gpu_utils/gpu_utils.h"
62 #include "gromacs/imd/imd.h"
63 #include "gromacs/listed_forces/bonded.h"
64 #include "gromacs/listed_forces/disre.h"
65 #include "gromacs/listed_forces/gpubonded.h"
66 #include "gromacs/listed_forces/listed_forces.h"
67 #include "gromacs/listed_forces/manage_threading.h"
68 #include "gromacs/listed_forces/orires.h"
69 #include "gromacs/math/arrayrefwithpadding.h"
70 #include "gromacs/math/functions.h"
71 #include "gromacs/math/units.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/math/vecdump.h"
74 #include "gromacs/mdlib/calcmu.h"
75 #include "gromacs/mdlib/calcvir.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/enerdata_utils.h"
78 #include "gromacs/mdlib/force.h"
79 #include "gromacs/mdlib/forcerec.h"
80 #include "gromacs/mdlib/gmx_omp_nthreads.h"
81 #include "gromacs/mdlib/ppforceworkload.h"
82 #include "gromacs/mdlib/qmmm.h"
83 #include "gromacs/mdlib/update.h"
84 #include "gromacs/mdtypes/commrec.h"
85 #include "gromacs/mdtypes/enerdata.h"
86 #include "gromacs/mdtypes/forceoutput.h"
87 #include "gromacs/mdtypes/iforceprovider.h"
88 #include "gromacs/mdtypes/inputrec.h"
89 #include "gromacs/mdtypes/md_enums.h"
90 #include "gromacs/mdtypes/state.h"
91 #include "gromacs/nbnxm/atomdata.h"
92 #include "gromacs/nbnxm/gpu_data_mgmt.h"
93 #include "gromacs/nbnxm/nbnxm.h"
94 #include "gromacs/pbcutil/ishift.h"
95 #include "gromacs/pbcutil/mshift.h"
96 #include "gromacs/pbcutil/pbc.h"
97 #include "gromacs/pulling/pull.h"
98 #include "gromacs/pulling/pull_rotation.h"
99 #include "gromacs/timing/cyclecounter.h"
100 #include "gromacs/timing/gpu_timing.h"
101 #include "gromacs/timing/wallcycle.h"
102 #include "gromacs/timing/wallcyclereporting.h"
103 #include "gromacs/timing/walltime_accounting.h"
104 #include "gromacs/topology/topology.h"
105 #include "gromacs/utility/arrayref.h"
106 #include "gromacs/utility/basedefinitions.h"
107 #include "gromacs/utility/cstringutil.h"
108 #include "gromacs/utility/exceptions.h"
109 #include "gromacs/utility/fatalerror.h"
110 #include "gromacs/utility/gmxassert.h"
111 #include "gromacs/utility/gmxmpi.h"
112 #include "gromacs/utility/logger.h"
113 #include "gromacs/utility/smalloc.h"
114 #include "gromacs/utility/strconvert.h"
115 #include "gromacs/utility/sysinfo.h"
116
117 // TODO: this environment variable allows us to verify before release
118 // that on less common architectures the total cost of polling is not larger than
119 // a blocking wait (so polling does not introduce overhead when the static
120 // PME-first ordering would suffice).
121 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
122
123
124 static void sum_forces(rvec f[], gmx::ArrayRef<const gmx::RVec> forceToAdd)
125 {
126     const int      end = forceToAdd.size();
127
128     int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
129 #pragma omp parallel for num_threads(nt) schedule(static)
130     for (int i = 0; i < end; i++)
131     {
132         rvec_inc(f[i], forceToAdd[i]);
133     }
134 }
135
136 static void calc_virial(int start, int homenr, const rvec x[], const rvec f[],
137                         tensor vir_part, const t_graph *graph, const matrix box,
138                         t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
139 {
140     /* The short-range virial from surrounding boxes */
141     calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
142     inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
143
144     /* Calculate partial virial, for local atoms only, based on short range.
145      * Total virial is computed in global_stat, called from do_md
146      */
147     f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
148     inc_nrnb(nrnb, eNR_VIRIAL, homenr);
149
150     if (debug)
151     {
152         pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
153     }
154 }
155
156 static void pull_potential_wrapper(const t_commrec *cr,
157                                    const t_inputrec *ir,
158                                    const matrix box, gmx::ArrayRef<const gmx::RVec> x,
159                                    gmx::ForceWithVirial *force,
160                                    const t_mdatoms *mdatoms,
161                                    gmx_enerdata_t *enerd,
162                                    pull_t *pull_work,
163                                    const real *lambda,
164                                    double t,
165                                    gmx_wallcycle_t wcycle)
166 {
167     t_pbc  pbc;
168     real   dvdl;
169
170     /* Calculate the center of mass forces, this requires communication,
171      * which is why pull_potential is called close to other communication.
172      */
173     wallcycle_start(wcycle, ewcPULLPOT);
174     set_pbc(&pbc, ir->ePBC, box);
175     dvdl                     = 0;
176     enerd->term[F_COM_PULL] +=
177         pull_potential(pull_work, mdatoms, &pbc,
178                        cr, t, lambda[efptRESTRAINT], as_rvec_array(x.data()), force, &dvdl);
179     enerd->dvdl_lin[efptRESTRAINT] += dvdl;
180     wallcycle_stop(wcycle, ewcPULLPOT);
181 }
182
183 static void pme_receive_force_ener(const t_commrec      *cr,
184                                    gmx::ForceWithVirial *forceWithVirial,
185                                    gmx_enerdata_t       *enerd,
186                                    gmx_wallcycle_t       wcycle)
187 {
188     real   e_q, e_lj, dvdl_q, dvdl_lj;
189     float  cycles_ppdpme, cycles_seppme;
190
191     cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
192     dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
193
194     /* In case of node-splitting, the PP nodes receive the long-range
195      * forces, virial and energy from the PME nodes here.
196      */
197     wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
198     dvdl_q  = 0;
199     dvdl_lj = 0;
200     gmx_pme_receive_f(cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
201                       &cycles_seppme);
202     enerd->term[F_COUL_RECIP] += e_q;
203     enerd->term[F_LJ_RECIP]   += e_lj;
204     enerd->dvdl_lin[efptCOUL] += dvdl_q;
205     enerd->dvdl_lin[efptVDW]  += dvdl_lj;
206
207     if (wcycle)
208     {
209         dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
210     }
211     wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
212 }
213
214 static void print_large_forces(FILE            *fp,
215                                const t_mdatoms *md,
216                                const t_commrec *cr,
217                                int64_t          step,
218                                real             forceTolerance,
219                                const rvec      *x,
220                                const rvec      *f)
221 {
222     real           force2Tolerance = gmx::square(forceTolerance);
223     gmx::index     numNonFinite    = 0;
224     for (int i = 0; i < md->homenr; i++)
225     {
226         real force2    = norm2(f[i]);
227         bool nonFinite = !std::isfinite(force2);
228         if (force2 >= force2Tolerance || nonFinite)
229         {
230             fprintf(fp, "step %" PRId64 " atom %6d  x %8.3f %8.3f %8.3f  force %12.5e\n",
231                     step,
232                     ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
233         }
234         if (nonFinite)
235         {
236             numNonFinite++;
237         }
238     }
239     if (numNonFinite > 0)
240     {
241         /* Note that with MPI this fatal call on one rank might interrupt
242          * the printing on other ranks. But we can only avoid that with
243          * an expensive MPI barrier that we would need at each step.
244          */
245         gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
246     }
247 }
248
249 static void post_process_forces(const t_commrec           *cr,
250                                 int64_t                    step,
251                                 t_nrnb                    *nrnb,
252                                 gmx_wallcycle_t            wcycle,
253                                 const gmx_localtop_t      *top,
254                                 const matrix               box,
255                                 const rvec                 x[],
256                                 rvec                       f[],
257                                 gmx::ForceWithVirial      *forceWithVirial,
258                                 tensor                     vir_force,
259                                 const t_mdatoms           *mdatoms,
260                                 const t_graph             *graph,
261                                 const t_forcerec          *fr,
262                                 const gmx_vsite_t         *vsite,
263                                 int                        flags)
264 {
265     if (fr->haveDirectVirialContributions)
266     {
267         rvec *fDirectVir = as_rvec_array(forceWithVirial->force_.data());
268
269         if (vsite)
270         {
271             /* Spread the mesh force on virtual sites to the other particles...
272              * This is parallellized. MPI communication is performed
273              * if the constructing atoms aren't local.
274              */
275             matrix virial = { { 0 } };
276             spread_vsite_f(vsite, x, fDirectVir, nullptr,
277                            (flags & GMX_FORCE_VIRIAL) != 0, virial,
278                            nrnb,
279                            &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
280             forceWithVirial->addVirialContribution(virial);
281         }
282
283         if (flags & GMX_FORCE_VIRIAL)
284         {
285             /* Now add the forces, this is local */
286             sum_forces(f, forceWithVirial->force_);
287
288             /* Add the direct virial contributions */
289             GMX_ASSERT(forceWithVirial->computeVirial_, "forceWithVirial should request virial computation when we request the virial");
290             m_add(vir_force, forceWithVirial->getVirial(), vir_force);
291
292             if (debug)
293             {
294                 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
295             }
296         }
297     }
298
299     if (fr->print_force >= 0)
300     {
301         print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
302     }
303 }
304
305 static void do_nb_verlet(t_forcerec                       *fr,
306                          const interaction_const_t        *ic,
307                          gmx_enerdata_t                   *enerd,
308                          const int                         flags,
309                          const Nbnxm::InteractionLocality  ilocality,
310                          const int                         clearF,
311                          const int64_t                     step,
312                          t_nrnb                           *nrnb,
313                          gmx_wallcycle_t                   wcycle)
314 {
315     if (!(flags & GMX_FORCE_NONBONDED))
316     {
317         /* skip non-bonded calculation */
318         return;
319     }
320
321     nonbonded_verlet_t *nbv  = fr->nbv.get();
322
323     /* GPU kernel launch overhead is already timed separately */
324     if (fr->cutoff_scheme != ecutsVERLET)
325     {
326         gmx_incons("Invalid cut-off scheme passed!");
327     }
328
329     if (!nbv->useGpu())
330     {
331         /* When dynamic pair-list  pruning is requested, we need to prune
332          * at nstlistPrune steps.
333          */
334         if (nbv->isDynamicPruningStepCpu(step))
335         {
336             /* Prune the pair-list beyond fr->ic->rlistPrune using
337              * the current coordinates of the atoms.
338              */
339             wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
340             nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
341             wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
342         }
343
344         wallcycle_sub_start(wcycle, ewcsNONBONDED);
345     }
346
347     nbv->dispatchNonbondedKernel(ilocality, *ic, flags, clearF, fr, enerd, nrnb);
348
349     if (!nbv->useGpu())
350     {
351         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
352     }
353 }
354
355 static inline void clear_rvecs_omp(int n, rvec v[])
356 {
357     int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, n);
358
359     /* Note that we would like to avoid this conditional by putting it
360      * into the omp pragma instead, but then we still take the full
361      * omp parallel for overhead (at least with gcc5).
362      */
363     if (nth == 1)
364     {
365         for (int i = 0; i < n; i++)
366         {
367             clear_rvec(v[i]);
368         }
369     }
370     else
371     {
372 #pragma omp parallel for num_threads(nth) schedule(static)
373         for (int i = 0; i < n; i++)
374         {
375             clear_rvec(v[i]);
376         }
377     }
378 }
379
380 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
381  *
382  * \param groupOptions  Group options, containing T-coupling options
383  */
384 static real averageKineticEnergyEstimate(const t_grpopts &groupOptions)
385 {
386     real nrdfCoupled   = 0;
387     real nrdfUncoupled = 0;
388     real kineticEnergy = 0;
389     for (int g = 0; g < groupOptions.ngtc; g++)
390     {
391         if (groupOptions.tau_t[g] >= 0)
392         {
393             nrdfCoupled   += groupOptions.nrdf[g];
394             kineticEnergy += groupOptions.nrdf[g]*0.5*groupOptions.ref_t[g]*BOLTZ;
395         }
396         else
397         {
398             nrdfUncoupled += groupOptions.nrdf[g];
399         }
400     }
401
402     /* This conditional with > also catches nrdf=0 */
403     if (nrdfCoupled > nrdfUncoupled)
404     {
405         return kineticEnergy*(nrdfCoupled + nrdfUncoupled)/nrdfCoupled;
406     }
407     else
408     {
409         return 0;
410     }
411 }
412
413 /*! \brief This routine checks that the potential energy is finite.
414  *
415  * Always checks that the potential energy is finite. If step equals
416  * inputrec.init_step also checks that the magnitude of the potential energy
417  * is reasonable. Terminates with a fatal error when a check fails.
418  * Note that passing this check does not guarantee finite forces,
419  * since those use slightly different arithmetics. But in most cases
420  * there is just a narrow coordinate range where forces are not finite
421  * and energies are finite.
422  *
423  * \param[in] step      The step number, used for checking and printing
424  * \param[in] enerd     The energy data; the non-bonded group energies need to be added to enerd.term[F_EPOT] before calling this routine
425  * \param[in] inputrec  The input record
426  */
427 static void checkPotentialEnergyValidity(int64_t               step,
428                                          const gmx_enerdata_t &enerd,
429                                          const t_inputrec     &inputrec)
430 {
431     /* Threshold valid for comparing absolute potential energy against
432      * the kinetic energy. Normally one should not consider absolute
433      * potential energy values, but with a factor of one million
434      * we should never get false positives.
435      */
436     constexpr real c_thresholdFactor = 1e6;
437
438     bool           energyIsNotFinite    = !std::isfinite(enerd.term[F_EPOT]);
439     real           averageKineticEnergy = 0;
440     /* We only check for large potential energy at the initial step,
441      * because that is by far the most likely step for this too occur
442      * and because computing the average kinetic energy is not free.
443      * Note: nstcalcenergy >> 1 often does not allow to catch large energies
444      * before they become NaN.
445      */
446     if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
447     {
448         averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
449     }
450
451     if (energyIsNotFinite || (averageKineticEnergy > 0 &&
452                               enerd.term[F_EPOT] > c_thresholdFactor*averageKineticEnergy))
453     {
454         gmx_fatal(FARGS, "Step %" PRId64 ": The total potential energy is %g, which is %s. The LJ and electrostatic contributions to the energy are %g and %g, respectively. A %s potential energy can be caused by overlapping interactions in bonded interactions or very large%s coordinate values. Usually this is caused by a badly- or non-equilibrated initial configuration, incorrect interactions or parameters in the topology.",
455                   step,
456                   enerd.term[F_EPOT],
457                   energyIsNotFinite ? "not finite" : "extremely high",
458                   enerd.term[F_LJ],
459                   enerd.term[F_COUL_SR],
460                   energyIsNotFinite ? "non-finite" : "very high",
461                   energyIsNotFinite ? " or Nan" : "");
462     }
463 }
464
465 /*! \brief Return true if there are special forces computed this step.
466  *
467  * The conditionals exactly correspond to those in computeSpecialForces().
468  */
469 static bool
470 haveSpecialForces(const t_inputrec              *inputrec,
471                   ForceProviders                *forceProviders,
472                   const pull_t                  *pull_work,
473                   int                            forceFlags,
474                   const gmx_edsam               *ed)
475 {
476     const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
477
478     return
479         ((computeForces && forceProviders->hasForceProvider()) ||         // forceProviders
480          (inputrec->bPull && pull_have_potential(pull_work)) ||           // pull
481          inputrec->bRot ||                                                // enforced rotation
482          (ed != nullptr) ||                                               // flooding
483          (inputrec->bIMD && computeForces));                              // IMD
484 }
485
486 /*! \brief Compute forces and/or energies for special algorithms
487  *
488  * The intention is to collect all calls to algorithms that compute
489  * forces on local atoms only and that do not contribute to the local
490  * virial sum (but add their virial contribution separately).
491  * Eventually these should likely all become ForceProviders.
492  * Within this function the intention is to have algorithms that do
493  * global communication at the end, so global barriers within the MD loop
494  * are as close together as possible.
495  *
496  * \param[in]     fplog            The log file
497  * \param[in]     cr               The communication record
498  * \param[in]     inputrec         The input record
499  * \param[in]     awh              The Awh module (nullptr if none in use).
500  * \param[in]     enforcedRotation Enforced rotation module.
501  * \param[in]     imdSession       The IMD session
502  * \param[in]     pull_work        The pull work structure.
503  * \param[in]     step             The current MD step
504  * \param[in]     t                The current time
505  * \param[in,out] wcycle           Wallcycle accounting struct
506  * \param[in,out] forceProviders   Pointer to a list of force providers
507  * \param[in]     box              The unit cell
508  * \param[in]     x                The coordinates
509  * \param[in]     mdatoms          Per atom properties
510  * \param[in]     lambda           Array of free-energy lambda values
511  * \param[in]     forceFlags       Flags that tell whether we should compute forces/energies/virial
512  * \param[in,out] forceWithVirial  Force and virial buffers
513  * \param[in,out] enerd            Energy buffer
514  * \param[in,out] ed               Essential dynamics pointer
515  * \param[in]     bNS              Tells if we did neighbor searching this step, used for ED sampling
516  *
517  * \todo Remove bNS, which is used incorrectly.
518  * \todo Convert all other algorithms called here to ForceProviders.
519  */
520 static void
521 computeSpecialForces(FILE                          *fplog,
522                      const t_commrec               *cr,
523                      const t_inputrec              *inputrec,
524                      gmx::Awh                      *awh,
525                      gmx_enfrot                    *enforcedRotation,
526                      gmx::ImdSession               *imdSession,
527                      pull_t                        *pull_work,
528                      int64_t                        step,
529                      double                         t,
530                      gmx_wallcycle_t                wcycle,
531                      ForceProviders                *forceProviders,
532                      matrix                         box,
533                      gmx::ArrayRef<const gmx::RVec> x,
534                      const t_mdatoms               *mdatoms,
535                      real                          *lambda,
536                      int                            forceFlags,
537                      gmx::ForceWithVirial          *forceWithVirial,
538                      gmx_enerdata_t                *enerd,
539                      gmx_edsam                     *ed,
540                      gmx_bool                       bNS)
541 {
542     const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
543
544     /* NOTE: Currently all ForceProviders only provide forces.
545      *       When they also provide energies, remove this conditional.
546      */
547     if (computeForces)
548     {
549         gmx::ForceProviderInput  forceProviderInput(x, *mdatoms, t, box, *cr);
550         gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
551
552         /* Collect forces from modules */
553         forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
554     }
555
556     if (inputrec->bPull && pull_have_potential(pull_work))
557     {
558         pull_potential_wrapper(cr, inputrec, box, x,
559                                forceWithVirial,
560                                mdatoms, enerd, pull_work, lambda, t,
561                                wcycle);
562
563         if (awh)
564         {
565             enerd->term[F_COM_PULL] +=
566                 awh->applyBiasForcesAndUpdateBias(inputrec->ePBC, *mdatoms, box,
567                                                   forceWithVirial,
568                                                   t, step, wcycle, fplog);
569         }
570     }
571
572     rvec *f = as_rvec_array(forceWithVirial->force_.data());
573
574     /* Add the forces from enforced rotation potentials (if any) */
575     if (inputrec->bRot)
576     {
577         wallcycle_start(wcycle, ewcROTadd);
578         enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
579         wallcycle_stop(wcycle, ewcROTadd);
580     }
581
582     if (ed)
583     {
584         /* Note that since init_edsam() is called after the initialization
585          * of forcerec, edsam doesn't request the noVirSum force buffer.
586          * Thus if no other algorithm (e.g. PME) requires it, the forces
587          * here will contribute to the virial.
588          */
589         do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, bNS);
590     }
591
592     /* Add forces from interactive molecular dynamics (IMD), if any */
593     if (inputrec->bIMD && computeForces)
594     {
595         imdSession->applyForces(f);
596     }
597 }
598
599 /*! \brief Launch the prepare_step and spread stages of PME GPU.
600  *
601  * \param[in]  pmedata       The PME structure
602  * \param[in]  box           The box matrix
603  * \param[in]  x             Coordinate array
604  * \param[in]  flags         Force flags
605  * \param[in]  pmeFlags      PME flags
606  * \param[in]  wcycle        The wallcycle structure
607  */
608 static inline void launchPmeGpuSpread(gmx_pme_t      *pmedata,
609                                       matrix          box,
610                                       rvec            x[],
611                                       int             flags,
612                                       int             pmeFlags,
613                                       gmx_wallcycle_t wcycle)
614 {
615     pme_gpu_prepare_computation(pmedata, (flags & GMX_FORCE_DYNAMICBOX) != 0, box, wcycle, pmeFlags);
616     pme_gpu_launch_spread(pmedata, x, wcycle);
617 }
618
619 /*! \brief Launch the FFT and gather stages of PME GPU
620  *
621  * This function only implements setting the output forces (no accumulation).
622  *
623  * \param[in]  pmedata        The PME structure
624  * \param[in]  wcycle         The wallcycle structure
625  */
626 static void launchPmeGpuFftAndGather(gmx_pme_t        *pmedata,
627                                      gmx_wallcycle_t   wcycle)
628 {
629     pme_gpu_launch_complex_transforms(pmedata, wcycle);
630     pme_gpu_launch_gather(pmedata, wcycle, PmeForceOutputHandling::Set);
631 }
632
633 /*! \brief
634  *  Polling wait for either of the PME or nonbonded GPU tasks.
635  *
636  * Instead of a static order in waiting for GPU tasks, this function
637  * polls checking which of the two tasks completes first, and does the
638  * associated force buffer reduction overlapped with the other task.
639  * By doing that, unlike static scheduling order, it can always overlap
640  * one of the reductions, regardless of the GPU task completion order.
641  *
642  * \param[in]     nbv              Nonbonded verlet structure
643  * \param[in,out] pmedata          PME module data
644  * \param[in,out] force            Force array to reduce task outputs into.
645  * \param[in,out] forceWithVirial  Force and virial buffers
646  * \param[in,out] fshift           Shift force output vector results are reduced into
647  * \param[in,out] enerd            Energy data structure results are reduced into
648  * \param[in]     flags            Force flags
649  * \param[in]     pmeFlags         PME flags
650  * \param[in]     haveOtherWork    Tells whether there is other work than non-bonded in the stream(s)
651  * \param[in]     wcycle           The wallcycle structure
652  */
653 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t                  *nbv,
654                                         gmx_pme_t                           *pmedata,
655                                         gmx::ArrayRefWithPadding<gmx::RVec> *force,
656                                         gmx::ForceWithVirial                *forceWithVirial,
657                                         rvec                                 fshift[],
658                                         gmx_enerdata_t                      *enerd,
659                                         int                                  flags,
660                                         int                                  pmeFlags,
661                                         bool                                 haveOtherWork,
662                                         gmx_wallcycle_t                      wcycle)
663 {
664     bool isPmeGpuDone = false;
665     bool isNbGpuDone  = false;
666
667
668     gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
669
670     while (!isPmeGpuDone || !isNbGpuDone)
671     {
672         if (!isPmeGpuDone)
673         {
674             GpuTaskCompletion completionType = (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
675             isPmeGpuDone = pme_gpu_try_finish_task(pmedata, pmeFlags, wcycle, forceWithVirial, enerd, completionType);
676         }
677
678         if (!isNbGpuDone)
679         {
680             GpuTaskCompletion completionType = (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
681             wallcycle_start_nocount(wcycle, ewcWAIT_GPU_NB_L);
682             isNbGpuDone = Nbnxm::gpu_try_finish_task(nbv->gpu_nbv,
683                                                      flags,
684                                                      Nbnxm::AtomLocality::Local,
685                                                      haveOtherWork,
686                                                      enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
687                                                      fshift, completionType);
688             wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
689             // To get the call count right, when the task finished we
690             // issue a start/stop.
691             // TODO: move the ewcWAIT_GPU_NB_L cycle counting into nbnxn_gpu_try_finish_task()
692             // and ewcNB_XF_BUF_OPS counting into nbnxn_atomdata_add_nbat_f_to_f().
693             if (isNbGpuDone)
694             {
695                 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
696                 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
697
698                 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
699                                               as_rvec_array(force->unpaddedArrayRef().data()), wcycle);
700             }
701         }
702     }
703 }
704
705 /*! \brief Hack structure with force ouput buffers for do_force for the home atoms for this domain */
706 struct ForceOutputs
707 {
708     //! Constructor
709     ForceOutputs(rvec *f, gmx::ForceWithVirial const forceWithVirial) :
710         f(f),
711         forceWithVirial(forceWithVirial) {}
712
713     //! Force output buffer used by legacy modules (without SIMD padding)
714     rvec                 *const f;
715     //! Force with direct virial contribution (if there are any; without SIMD padding)
716     gmx::ForceWithVirial        forceWithVirial;
717 };
718
719 /*! \brief Set up the different force buffers; also does clearing.
720  *
721  * \param[in] fr        force record pointer
722  * \param[in] pull_work The pull work object.
723  * \param[in] inputrec  input record
724  * \param[in] force     force array
725  * \param[in] bDoForces True if force are computed this step
726  * \param[in] doVirial  True if virial is computed this step
727  * \param[out] wcycle   wallcycle recording structure
728  *
729  * \returns             Cleared force output structure
730  */
731 static ForceOutputs
732 setupForceOutputs(const t_forcerec                    *fr,
733                   pull_t                              *pull_work,
734                   const t_inputrec                    &inputrec,
735                   gmx::ArrayRefWithPadding<gmx::RVec>  force,
736                   const bool                           bDoForces,
737                   const bool                           doVirial,
738                   gmx_wallcycle_t                      wcycle)
739 {
740     wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
741
742     /* Temporary solution until all routines take PaddedRVecVector */
743     rvec *const f = as_rvec_array(force.unpaddedArrayRef().data());
744     if (bDoForces)
745     {
746         /* Clear the short- and long-range forces */
747         clear_rvecs_omp(fr->natoms_force_constr, f);
748     }
749
750     /* If we need to compute the virial, we might need a separate
751      * force buffer for algorithms for which the virial is calculated
752      * directly, such as PME. Otherwise, forceWithVirial uses the
753      * the same force (f in legacy calls) buffer as other algorithms.
754      */
755     const bool useSeparateForceWithVirialBuffer = (bDoForces && (doVirial && fr->haveDirectVirialContributions));
756
757
758     /* forceWithVirial uses the local atom range only */
759     gmx::ForceWithVirial forceWithVirial (useSeparateForceWithVirialBuffer ?
760                                           *fr->forceBufferForDirectVirialContributions : force.unpaddedArrayRef(),
761                                           doVirial);
762
763     if (useSeparateForceWithVirialBuffer)
764     {
765         /* TODO: update comment
766          * We only compute forces on local atoms. Note that vsites can
767          * spread to non-local atoms, but that part of the buffer is
768          * cleared separately in the vsite spreading code.
769          */
770         clear_rvecs_omp(forceWithVirial.force_.size(), as_rvec_array(forceWithVirial.force_.data()));
771     }
772
773     if (inputrec.bPull && pull_have_constraint(pull_work))
774     {
775         clear_pull_forces(pull_work);
776     }
777
778     wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
779
780     return ForceOutputs(f, forceWithVirial);
781 }
782
783
784 /*! \brief Set up flags that indicate what type of work is there to compute.
785  *
786  * Currently we only update it at search steps,
787  * but some properties may change more frequently (e.g. virial/non-virial step),
788  * so when including those either the frequency of update (per-step) or the scope
789  * of a flag will change (i.e. a set of flags for nstlist steps).
790  *
791  */
792 static void
793 setupForceWorkload(gmx::PpForceWorkload *forceWork,
794                    const t_inputrec     *inputrec,
795                    const t_forcerec     *fr,
796                    const pull_t         *pull_work,
797                    const gmx_edsam      *ed,
798                    const t_idef         &idef,
799                    const t_fcdata       *fcd,
800                    const int             forceFlags
801                    )
802 {
803     forceWork->haveSpecialForces      = haveSpecialForces(inputrec, fr->forceProviders, pull_work, forceFlags, ed);
804     forceWork->haveCpuBondedWork      = haveCpuBondeds(*fr);
805     forceWork->haveGpuBondedWork      = ((fr->gpuBonded != nullptr) && fr->gpuBonded->haveInteractions());
806     forceWork->haveRestraintsWork     = havePositionRestraints(idef, *fcd);
807     forceWork->haveCpuListedForceWork = haveCpuListedForces(*fr, idef, *fcd);
808 }
809
810 void do_force(FILE                                     *fplog,
811               const t_commrec                          *cr,
812               const gmx_multisim_t                     *ms,
813               const t_inputrec                         *inputrec,
814               gmx::Awh                                 *awh,
815               gmx_enfrot                               *enforcedRotation,
816               gmx::ImdSession                          *imdSession,
817               pull_t                                   *pull_work,
818               int64_t                                   step,
819               t_nrnb                                   *nrnb,
820               gmx_wallcycle_t                           wcycle,
821               const gmx_localtop_t                     *top,
822               matrix                                    box,
823               gmx::ArrayRefWithPadding<gmx::RVec>       x,     //NOLINT(performance-unnecessary-value-param)
824               history_t                                *hist,
825               gmx::ArrayRefWithPadding<gmx::RVec>       force, //NOLINT(performance-unnecessary-value-param)
826               tensor                                    vir_force,
827               const t_mdatoms                          *mdatoms,
828               gmx_enerdata_t                           *enerd,
829               t_fcdata                                 *fcd,
830               gmx::ArrayRef<real>                       lambda,
831               t_graph                                  *graph,
832               t_forcerec                               *fr,
833               gmx::PpForceWorkload                     *ppForceWorkload,
834               const gmx_vsite_t                        *vsite,
835               rvec                                      mu_tot,
836               double                                    t,
837               gmx_edsam                                *ed,
838               int                                       flags,
839               const DDBalanceRegionHandler             &ddBalanceRegionHandler)
840 {
841     int                  i, j;
842     double               mu[2*DIM];
843     gmx_bool             bStateChanged, bNS, bFillGrid, bCalcCGCM;
844     gmx_bool             bDoForces, bUseGPU, bUseOrEmulGPU;
845     nonbonded_verlet_t  *nbv = fr->nbv.get();
846     interaction_const_t *ic  = fr->ic;
847
848     /* modify force flag if not doing nonbonded */
849     if (!fr->bNonbonded)
850     {
851         flags &= ~GMX_FORCE_NONBONDED;
852     }
853     bStateChanged = ((flags & GMX_FORCE_STATECHANGED) != 0);
854     bNS           = ((flags & GMX_FORCE_NS) != 0);
855     bFillGrid     = (bNS && bStateChanged);
856     bCalcCGCM     = (bFillGrid && !DOMAINDECOMP(cr));
857     bDoForces     = ((flags & GMX_FORCE_FORCES) != 0);
858     bUseGPU       = fr->nbv->useGpu();
859     bUseOrEmulGPU = bUseGPU || fr->nbv->emulateGpu();
860
861     const auto pmeRunMode = fr->pmedata ? pme_run_mode(fr->pmedata) : PmeRunMode::CPU;
862     // TODO slim this conditional down - inputrec and duty checks should mean the same in proper code!
863     const bool useGpuPme  = EEL_PME(fr->ic->eeltype) && thisRankHasDuty(cr, DUTY_PME) &&
864         ((pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed));
865     const int  pmeFlags = GMX_PME_SPREAD | GMX_PME_SOLVE |
866         ((flags & GMX_FORCE_VIRIAL) ? GMX_PME_CALC_ENER_VIR : 0) |
867         ((flags & GMX_FORCE_ENERGY) ? GMX_PME_CALC_ENER_VIR : 0) |
868         ((flags & GMX_FORCE_FORCES) ? GMX_PME_CALC_F : 0);
869
870     /* At a search step we need to start the first balancing region
871      * somewhere early inside the step after communication during domain
872      * decomposition (and not during the previous step as usual).
873      */
874     if (bNS)
875     {
876         ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
877     }
878
879     const int start  = 0;
880     const int homenr = mdatoms->homenr;
881
882     clear_mat(vir_force);
883
884     if (bStateChanged)
885     {
886         update_forcerec(fr, box);
887
888         if (inputrecNeedMutot(inputrec))
889         {
890             /* Calculate total (local) dipole moment in a temporary common array.
891              * This makes it possible to sum them over nodes faster.
892              */
893             calc_mu(start, homenr,
894                     x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
895                     mu, mu+DIM);
896         }
897     }
898
899     if (fr->ePBC != epbcNONE)
900     {
901         /* Compute shift vectors every step,
902          * because of pressure coupling or box deformation!
903          */
904         if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
905         {
906             calc_shifts(box, fr->shift_vec);
907         }
908
909         if (bCalcCGCM)
910         {
911             put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr), gmx_omp_nthreads_get(emntDefault));
912             inc_nrnb(nrnb, eNR_SHIFTX, homenr);
913         }
914         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
915         {
916             unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
917         }
918     }
919
920     nbnxn_atomdata_copy_shiftvec((flags & GMX_FORCE_DYNAMICBOX) != 0,
921                                  fr->shift_vec, nbv->nbat.get());
922
923 #if GMX_MPI
924     if (!thisRankHasDuty(cr, DUTY_PME))
925     {
926         /* Send particle coordinates to the pme nodes.
927          * Since this is only implemented for domain decomposition
928          * and domain decomposition does not use the graph,
929          * we do not need to worry about shifting.
930          */
931         gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
932                                  lambda[efptCOUL], lambda[efptVDW],
933                                  (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)) != 0,
934                                  step, wcycle);
935     }
936 #endif /* GMX_MPI */
937
938     if (useGpuPme)
939     {
940         launchPmeGpuSpread(fr->pmedata, box, as_rvec_array(x.unpaddedArrayRef().data()), flags, pmeFlags, wcycle);
941     }
942
943     /* do gridding for pair search */
944     if (bNS)
945     {
946         if (graph && bStateChanged)
947         {
948             /* Calculate intramolecular shift vectors to make molecules whole */
949             mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
950         }
951
952         // TODO
953         // - vzero is constant, do we need to pass it?
954         // - box_diag should be passed directly to nbnxn_put_on_grid
955         //
956         rvec vzero;
957         clear_rvec(vzero);
958
959         rvec box_diag;
960         box_diag[XX] = box[XX][XX];
961         box_diag[YY] = box[YY][YY];
962         box_diag[ZZ] = box[ZZ][ZZ];
963
964         wallcycle_start(wcycle, ewcNS);
965         if (!DOMAINDECOMP(cr))
966         {
967             wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
968             nbnxn_put_on_grid(nbv, box,
969                               0, vzero, box_diag,
970                               nullptr, 0, mdatoms->homenr, -1,
971                               fr->cginfo, x.unpaddedArrayRef(),
972                               0, nullptr);
973             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
974         }
975         else
976         {
977             wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
978             nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd),
979                                        fr->cginfo, x.unpaddedArrayRef());
980             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
981         }
982
983         nbv->setAtomProperties(*mdatoms, *fr->cginfo);
984
985         wallcycle_stop(wcycle, ewcNS);
986
987         /* initialize the GPU nbnxm atom data and bonded data structures */
988         if (bUseGPU)
989         {
990             wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
991
992             wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
993             Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
994             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
995
996             if (fr->gpuBonded)
997             {
998                 /* Now we put all atoms on the grid, we can assign bonded
999                  * interactions to the GPU, where the grid order is
1000                  * needed. Also the xq, f and fshift device buffers have
1001                  * been reallocated if needed, so the bonded code can
1002                  * learn about them. */
1003                 // TODO the xq, f, and fshift buffers are now shared
1004                 // resources, so they should be maintained by a
1005                 // higher-level object than the nb module.
1006                 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1007                                                                       top->idef,
1008                                                                       Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1009                                                                       Nbnxm::gpu_get_f(nbv->gpu_nbv),
1010                                                                       Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1011             }
1012             wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1013         }
1014
1015         // Need to run after the GPU-offload bonded interaction lists
1016         // are set up to be able to determine whether there is bonded work.
1017         setupForceWorkload(ppForceWorkload,
1018                            inputrec,
1019                            fr,
1020                            pull_work,
1021                            ed,
1022                            top->idef,
1023                            fcd,
1024                            flags);
1025     }
1026
1027     /* do local pair search */
1028     if (bNS)
1029     {
1030         // TODO: fuse this branch with the above bNS block
1031         wallcycle_start_nocount(wcycle, ewcNS);
1032         wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1033         /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1034         nbv->constructPairlist(Nbnxm::InteractionLocality::Local,
1035                                &top->excls, step, nrnb);
1036         wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1037         wallcycle_stop(wcycle, ewcNS);
1038     }
1039     else
1040     {
1041         nbv->setCoordinates(Nbnxm::AtomLocality::Local, false,
1042                             x.unpaddedArrayRef(), wcycle);
1043     }
1044
1045     if (bUseGPU)
1046     {
1047         ddBalanceRegionHandler.openBeforeForceComputationGpu();
1048
1049         wallcycle_start(wcycle, ewcLAUNCH_GPU);
1050
1051         wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1052         Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1053         Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1054                                   Nbnxm::AtomLocality::Local,
1055                                   ppForceWorkload->haveGpuBondedWork);
1056         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1057
1058         // bonded work not split into separate local and non-local, so with DD
1059         // we can only launch the kernel after non-local coordinates have been received.
1060         if (ppForceWorkload->haveGpuBondedWork && !havePPDomainDecomposition(cr))
1061         {
1062             wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1063             fr->gpuBonded->launchKernels(fr, flags, box);
1064             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1065         }
1066
1067         /* launch local nonbonded work on GPU */
1068         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1069         do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFNo,
1070                      step, nrnb, wcycle);
1071         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1072         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1073     }
1074
1075     if (useGpuPme)
1076     {
1077         // In PME GPU and mixed mode we launch FFT / gather after the
1078         // X copy/transform to allow overlap as well as after the GPU NB
1079         // launch to avoid FFT launch overhead hijacking the CPU and delaying
1080         // the nonbonded kernel.
1081         launchPmeGpuFftAndGather(fr->pmedata, wcycle);
1082     }
1083
1084     /* Communicate coordinates and sum dipole if necessary +
1085        do non-local pair search */
1086     if (havePPDomainDecomposition(cr))
1087     {
1088         if (bNS)
1089         {
1090             // TODO: fuse this branch with the above large bNS block
1091             wallcycle_start_nocount(wcycle, ewcNS);
1092             wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1093             /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1094             nbv->constructPairlist(Nbnxm::InteractionLocality::NonLocal,
1095                                    &top->excls, step, nrnb);
1096             wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1097             wallcycle_stop(wcycle, ewcNS);
1098         }
1099         else
1100         {
1101             dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1102
1103             nbv->setCoordinates(Nbnxm::AtomLocality::NonLocal, false,
1104                                 x.unpaddedArrayRef(), wcycle);
1105         }
1106
1107         if (bUseGPU)
1108         {
1109             wallcycle_start(wcycle, ewcLAUNCH_GPU);
1110
1111             /* launch non-local nonbonded tasks on GPU */
1112             wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1113             Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1114                                       Nbnxm::AtomLocality::NonLocal,
1115                                       ppForceWorkload->haveGpuBondedWork);
1116             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1117
1118             if (ppForceWorkload->haveGpuBondedWork)
1119             {
1120                 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1121                 fr->gpuBonded->launchKernels(fr, flags, box);
1122                 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1123             }
1124
1125             wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1126             do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1127                          step, nrnb, wcycle);
1128             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1129
1130             wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1131         }
1132     }
1133
1134     if (bUseGPU)
1135     {
1136         /* launch D2H copy-back F */
1137         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1138         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1139         if (havePPDomainDecomposition(cr))
1140         {
1141             Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1142                                       flags, Nbnxm::AtomLocality::NonLocal, ppForceWorkload->haveGpuBondedWork);
1143         }
1144         Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1145                                   flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork);
1146         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1147
1148         if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1149         {
1150             wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1151             fr->gpuBonded->launchEnergyTransfer();
1152             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1153         }
1154         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1155     }
1156
1157     if (bStateChanged && inputrecNeedMutot(inputrec))
1158     {
1159         if (PAR(cr))
1160         {
1161             gmx_sumd(2*DIM, mu, cr);
1162
1163             ddBalanceRegionHandler.reopenRegionCpu();
1164         }
1165
1166         for (i = 0; i < 2; i++)
1167         {
1168             for (j = 0; j < DIM; j++)
1169             {
1170                 fr->mu_tot[i][j] = mu[i*DIM + j];
1171             }
1172         }
1173     }
1174     if (fr->efep == efepNO)
1175     {
1176         copy_rvec(fr->mu_tot[0], mu_tot);
1177     }
1178     else
1179     {
1180         for (j = 0; j < DIM; j++)
1181         {
1182             mu_tot[j] =
1183                 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1184                 lambda[efptCOUL]*fr->mu_tot[1][j];
1185         }
1186     }
1187
1188     /* Reset energies */
1189     reset_enerdata(enerd);
1190     clear_rvecs(SHIFTS, fr->fshift);
1191
1192     if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1193     {
1194         wallcycle_start(wcycle, ewcPPDURINGPME);
1195         dd_force_flop_start(cr->dd, nrnb);
1196     }
1197
1198     if (inputrec->bRot)
1199     {
1200         wallcycle_start(wcycle, ewcROT);
1201         do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, bNS);
1202         wallcycle_stop(wcycle, ewcROT);
1203     }
1204
1205     /* Start the force cycle counter.
1206      * Note that a different counter is used for dynamic load balancing.
1207      */
1208     wallcycle_start(wcycle, ewcFORCE);
1209
1210     // set up and clear force outputs
1211     struct ForceOutputs forceOut = setupForceOutputs(fr, pull_work, *inputrec, force, bDoForces,
1212                                                      ((flags & GMX_FORCE_VIRIAL) != 0), wcycle);
1213
1214     /* We calculate the non-bonded forces, when done on the CPU, here.
1215      * We do this before calling do_force_lowlevel, because in that
1216      * function, the listed forces are calculated before PME, which
1217      * does communication.  With this order, non-bonded and listed
1218      * force calculation imbalance can be balanced out by the domain
1219      * decomposition load balancing.
1220      */
1221
1222     if (!bUseOrEmulGPU)
1223     {
1224         do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFYes,
1225                      step, nrnb, wcycle);
1226     }
1227
1228     if (fr->efep != efepNO)
1229     {
1230         /* Calculate the local and non-local free energy interactions here.
1231          * Happens here on the CPU both with and without GPU.
1232          */
1233         wallcycle_sub_start(wcycle, ewcsNONBONDED);
1234         nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::Local,
1235                                       fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1236                                       inputrec->fepvals, lambda.data(),
1237                                       enerd, flags, nrnb);
1238
1239         if (havePPDomainDecomposition(cr))
1240         {
1241             nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::NonLocal,
1242                                           fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1243                                           inputrec->fepvals, lambda.data(),
1244                                           enerd, flags, nrnb);
1245         }
1246         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
1247     }
1248
1249     if (!bUseOrEmulGPU)
1250     {
1251         if (havePPDomainDecomposition(cr))
1252         {
1253             do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1254                          step, nrnb, wcycle);
1255         }
1256
1257         /* Add all the non-bonded force to the normal force array.
1258          * This can be split into a local and a non-local part when overlapping
1259          * communication with calculation with domain decomposition.
1260          */
1261         wallcycle_stop(wcycle, ewcFORCE);
1262
1263         nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::All, forceOut.f, wcycle);
1264
1265         wallcycle_start_nocount(wcycle, ewcFORCE);
1266
1267         /* If there are multiple fshift output buffers we need to reduce them */
1268         if (flags & GMX_FORCE_VIRIAL)
1269         {
1270             /* This is not in a subcounter because it takes a
1271                negligible and constant-sized amount of time */
1272             nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->nbat.get(),
1273                                                      fr->fshift);
1274         }
1275     }
1276
1277     /* update QMMMrec, if necessary */
1278     if (fr->bQMMM)
1279     {
1280         update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1281     }
1282
1283     /* Compute the bonded and non-bonded energies and optionally forces */
1284     do_force_lowlevel(fr, inputrec, &(top->idef),
1285                       cr, ms, nrnb, wcycle, mdatoms,
1286                       as_rvec_array(x.unpaddedArrayRef().data()), hist, forceOut.f, &forceOut.forceWithVirial, enerd, fcd,
1287                       box, lambda.data(), graph, &(top->excls), fr->mu_tot,
1288                       flags,
1289                       ddBalanceRegionHandler);
1290
1291     wallcycle_stop(wcycle, ewcFORCE);
1292
1293     computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1294                          imdSession, pull_work, step, t, wcycle,
1295                          fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda.data(),
1296                          flags, &forceOut.forceWithVirial, enerd,
1297                          ed, bNS);
1298
1299     // Will store the amount of cycles spent waiting for the GPU that
1300     // will be later used in the DLB accounting.
1301     float cycles_wait_gpu = 0;
1302     if (bUseOrEmulGPU)
1303     {
1304         /* wait for non-local forces (or calculate in emulation mode) */
1305         if (havePPDomainDecomposition(cr))
1306         {
1307             if (bUseGPU)
1308             {
1309                 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1310                 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1311                                             flags, Nbnxm::AtomLocality::NonLocal,
1312                                             ppForceWorkload->haveGpuBondedWork,
1313                                             enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1314                                             fr->fshift);
1315                 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1316             }
1317             else
1318             {
1319                 wallcycle_start_nocount(wcycle, ewcFORCE);
1320                 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFYes,
1321                              step, nrnb, wcycle);
1322                 wallcycle_stop(wcycle, ewcFORCE);
1323             }
1324
1325             nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::NonLocal,
1326                                           forceOut.f, wcycle);
1327         }
1328     }
1329
1330     if (havePPDomainDecomposition(cr))
1331     {
1332         /* We are done with the CPU compute.
1333          * We will now communicate the non-local forces.
1334          * If we use a GPU this will overlap with GPU work, so in that case
1335          * we do not close the DD force balancing region here.
1336          */
1337         ddBalanceRegionHandler.closeAfterForceComputationCpu();
1338
1339         if (bDoForces)
1340         {
1341             dd_move_f(cr->dd, force.unpaddedArrayRef(), fr->fshift, wcycle);
1342         }
1343     }
1344
1345     // With both nonbonded and PME offloaded a GPU on the same rank, we use
1346     // an alternating wait/reduction scheme.
1347     bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPme && bUseGPU && !DOMAINDECOMP(cr));
1348     if (alternateGpuWait)
1349     {
1350         alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &force, &forceOut.forceWithVirial, fr->fshift, enerd,
1351                                     flags, pmeFlags, ppForceWorkload->haveGpuBondedWork, wcycle);
1352     }
1353
1354     if (!alternateGpuWait && useGpuPme)
1355     {
1356         pme_gpu_wait_and_reduce(fr->pmedata, pmeFlags, wcycle, &forceOut.forceWithVirial, enerd);
1357     }
1358
1359     /* Wait for local GPU NB outputs on the non-alternating wait path */
1360     if (!alternateGpuWait && bUseGPU)
1361     {
1362         /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1363          * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1364          * but even with a step of 0.1 ms the difference is less than 1%
1365          * of the step time.
1366          */
1367         const float gpuWaitApiOverheadMargin = 2e6f; /* cycles */
1368
1369         wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1370         Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1371                                     flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork,
1372                                     enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1373                                     fr->fshift);
1374         float cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1375
1376         if (ddBalanceRegionHandler.useBalancingRegion())
1377         {
1378             DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1379             if (bDoForces && cycles_tmp <= gpuWaitApiOverheadMargin)
1380             {
1381                 /* We measured few cycles, it could be that the kernel
1382                  * and transfer finished earlier and there was no actual
1383                  * wait time, only API call overhead.
1384                  * Then the actual time could be anywhere between 0 and
1385                  * cycles_wait_est. We will use half of cycles_wait_est.
1386                  */
1387                 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1388             }
1389             ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1390         }
1391     }
1392
1393     if (fr->nbv->emulateGpu())
1394     {
1395         // NOTE: emulation kernel is not included in the balancing region,
1396         // but emulation mode does not target performance anyway
1397         wallcycle_start_nocount(wcycle, ewcFORCE);
1398         do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local,
1399                      DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1400                      step, nrnb, wcycle);
1401         wallcycle_stop(wcycle, ewcFORCE);
1402     }
1403
1404     if (useGpuPme)
1405     {
1406         pme_gpu_reinit_computation(fr->pmedata, wcycle);
1407     }
1408
1409     if (bUseGPU)
1410     {
1411         /* now clear the GPU outputs while we finish the step on the CPU */
1412         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1413         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1414         Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, flags);
1415
1416         if (nbv->isDynamicPruningStepGpu(step))
1417         {
1418             nbv->dispatchPruneKernelGpu(step);
1419         }
1420         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1421         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1422     }
1423
1424     if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1425     {
1426         wallcycle_start(wcycle, ewcWAIT_GPU_BONDED);
1427         // in principle this should be included in the DD balancing region,
1428         // but generally it is infrequent so we'll omit it for the sake of
1429         // simpler code
1430         fr->gpuBonded->accumulateEnergyTerms(enerd);
1431         wallcycle_stop(wcycle, ewcWAIT_GPU_BONDED);
1432
1433         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1434         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1435         fr->gpuBonded->clearEnergies();
1436         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1437         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1438     }
1439
1440     /* Do the nonbonded GPU (or emulation) force buffer reduction
1441      * on the non-alternating path. */
1442     if (bUseOrEmulGPU && !alternateGpuWait)
1443     {
1444         nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
1445                                       forceOut.f, wcycle);
1446     }
1447     if (DOMAINDECOMP(cr))
1448     {
1449         dd_force_flop_stop(cr->dd, nrnb);
1450     }
1451
1452     if (bDoForces)
1453     {
1454         /* If we have NoVirSum forces, but we do not calculate the virial,
1455          * we sum fr->f_novirsum=forceOut.f later.
1456          */
1457         if (vsite && !(fr->haveDirectVirialContributions && !(flags & GMX_FORCE_VIRIAL)))
1458         {
1459             spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, fr->fshift, FALSE, nullptr, nrnb,
1460                            &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1461         }
1462
1463         if (flags & GMX_FORCE_VIRIAL)
1464         {
1465             /* Calculation of the virial must be done after vsites! */
1466             calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f,
1467                         vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1468         }
1469     }
1470
1471     if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1472     {
1473         /* In case of node-splitting, the PP nodes receive the long-range
1474          * forces, virial and energy from the PME nodes here.
1475          */
1476         pme_receive_force_ener(cr, &forceOut.forceWithVirial, enerd, wcycle);
1477     }
1478
1479     if (bDoForces)
1480     {
1481         post_process_forces(cr, step, nrnb, wcycle,
1482                             top, box, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, &forceOut.forceWithVirial,
1483                             vir_force, mdatoms, graph, fr, vsite,
1484                             flags);
1485     }
1486
1487     if (flags & GMX_FORCE_ENERGY)
1488     {
1489         /* Sum the potential energy terms from group contributions */
1490         sum_epot(&(enerd->grpp), enerd->term);
1491
1492         if (!EI_TPI(inputrec->eI))
1493         {
1494             checkPotentialEnergyValidity(step, *enerd, *inputrec);
1495         }
1496     }
1497
1498     /* In case we don't have constraints and are using GPUs, the next balancing
1499      * region starts here.
1500      * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1501      * virial calculation and COM pulling, is not thus not included in
1502      * the balance timing, which is ok as most tasks do communication.
1503      */
1504     ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);
1505 }