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