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