Refactor gmx_group_t to SimulationAtomGroups
[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/force.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/gmx_omp_nthreads.h"
80 #include "gromacs/mdlib/ppforceworkload.h"
81 #include "gromacs/mdlib/qmmm.h"
82 #include "gromacs/mdlib/update.h"
83 #include "gromacs/mdtypes/commrec.h"
84 #include "gromacs/mdtypes/enerdata.h"
85 #include "gromacs/mdtypes/forceoutput.h"
86 #include "gromacs/mdtypes/iforceprovider.h"
87 #include "gromacs/mdtypes/inputrec.h"
88 #include "gromacs/mdtypes/md_enums.h"
89 #include "gromacs/mdtypes/state.h"
90 #include "gromacs/nbnxm/atomdata.h"
91 #include "gromacs/nbnxm/gpu_data_mgmt.h"
92 #include "gromacs/nbnxm/nbnxm.h"
93 #include "gromacs/pbcutil/ishift.h"
94 #include "gromacs/pbcutil/mshift.h"
95 #include "gromacs/pbcutil/pbc.h"
96 #include "gromacs/pulling/pull.h"
97 #include "gromacs/pulling/pull_rotation.h"
98 #include "gromacs/timing/cyclecounter.h"
99 #include "gromacs/timing/gpu_timing.h"
100 #include "gromacs/timing/wallcycle.h"
101 #include "gromacs/timing/wallcyclereporting.h"
102 #include "gromacs/timing/walltime_accounting.h"
103 #include "gromacs/topology/topology.h"
104 #include "gromacs/utility/arrayref.h"
105 #include "gromacs/utility/basedefinitions.h"
106 #include "gromacs/utility/cstringutil.h"
107 #include "gromacs/utility/exceptions.h"
108 #include "gromacs/utility/fatalerror.h"
109 #include "gromacs/utility/gmxassert.h"
110 #include "gromacs/utility/gmxmpi.h"
111 #include "gromacs/utility/logger.h"
112 #include "gromacs/utility/smalloc.h"
113 #include "gromacs/utility/strconvert.h"
114 #include "gromacs/utility/sysinfo.h"
115
116 // TODO: this environment variable allows us to verify before release
117 // that on less common architectures the total cost of polling is not larger than
118 // a blocking wait (so polling does not introduce overhead when the static
119 // PME-first ordering would suffice).
120 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
121
122
123 static void sum_forces(rvec f[], gmx::ArrayRef<const gmx::RVec> forceToAdd)
124 {
125     const int      end = forceToAdd.size();
126
127     int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
128 #pragma omp parallel for num_threads(nt) schedule(static)
129     for (int i = 0; i < end; i++)
130     {
131         rvec_inc(f[i], forceToAdd[i]);
132     }
133 }
134
135 static void calc_virial(int start, int homenr, const rvec x[], const rvec f[],
136                         tensor vir_part, const t_graph *graph, const matrix box,
137                         t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
138 {
139     /* The short-range virial from surrounding boxes */
140     calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
141     inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
142
143     /* Calculate partial virial, for local atoms only, based on short range.
144      * Total virial is computed in global_stat, called from do_md
145      */
146     f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
147     inc_nrnb(nrnb, eNR_VIRIAL, homenr);
148
149     if (debug)
150     {
151         pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
152     }
153 }
154
155 static void pull_potential_wrapper(const t_commrec *cr,
156                                    const t_inputrec *ir,
157                                    const matrix box, gmx::ArrayRef<const gmx::RVec> x,
158                                    gmx::ForceWithVirial *force,
159                                    const t_mdatoms *mdatoms,
160                                    gmx_enerdata_t *enerd,
161                                    const real *lambda,
162                                    double t,
163                                    gmx_wallcycle_t wcycle)
164 {
165     t_pbc  pbc;
166     real   dvdl;
167
168     /* Calculate the center of mass forces, this requires communication,
169      * which is why pull_potential is called close to other communication.
170      */
171     wallcycle_start(wcycle, ewcPULLPOT);
172     set_pbc(&pbc, ir->ePBC, box);
173     dvdl                     = 0;
174     enerd->term[F_COM_PULL] +=
175         pull_potential(ir->pull_work, mdatoms, &pbc,
176                        cr, t, lambda[efptRESTRAINT], as_rvec_array(x.data()), force, &dvdl);
177     enerd->dvdl_lin[efptRESTRAINT] += dvdl;
178     wallcycle_stop(wcycle, ewcPULLPOT);
179 }
180
181 static void pme_receive_force_ener(const t_commrec      *cr,
182                                    gmx::ForceWithVirial *forceWithVirial,
183                                    gmx_enerdata_t       *enerd,
184                                    gmx_wallcycle_t       wcycle)
185 {
186     real   e_q, e_lj, dvdl_q, dvdl_lj;
187     float  cycles_ppdpme, cycles_seppme;
188
189     cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
190     dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
191
192     /* In case of node-splitting, the PP nodes receive the long-range
193      * forces, virial and energy from the PME nodes here.
194      */
195     wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
196     dvdl_q  = 0;
197     dvdl_lj = 0;
198     gmx_pme_receive_f(cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
199                       &cycles_seppme);
200     enerd->term[F_COUL_RECIP] += e_q;
201     enerd->term[F_LJ_RECIP]   += e_lj;
202     enerd->dvdl_lin[efptCOUL] += dvdl_q;
203     enerd->dvdl_lin[efptVDW]  += dvdl_lj;
204
205     if (wcycle)
206     {
207         dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
208     }
209     wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
210 }
211
212 static void print_large_forces(FILE            *fp,
213                                const t_mdatoms *md,
214                                const t_commrec *cr,
215                                int64_t          step,
216                                real             forceTolerance,
217                                const rvec      *x,
218                                const rvec      *f)
219 {
220     real           force2Tolerance = gmx::square(forceTolerance);
221     gmx::index     numNonFinite    = 0;
222     for (int i = 0; i < md->homenr; i++)
223     {
224         real force2    = norm2(f[i]);
225         bool nonFinite = !std::isfinite(force2);
226         if (force2 >= force2Tolerance || nonFinite)
227         {
228             fprintf(fp, "step %" PRId64 " atom %6d  x %8.3f %8.3f %8.3f  force %12.5e\n",
229                     step,
230                     ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
231         }
232         if (nonFinite)
233         {
234             numNonFinite++;
235         }
236     }
237     if (numNonFinite > 0)
238     {
239         /* Note that with MPI this fatal call on one rank might interrupt
240          * the printing on other ranks. But we can only avoid that with
241          * an expensive MPI barrier that we would need at each step.
242          */
243         gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
244     }
245 }
246
247 static void post_process_forces(const t_commrec           *cr,
248                                 int64_t                    step,
249                                 t_nrnb                    *nrnb,
250                                 gmx_wallcycle_t            wcycle,
251                                 const gmx_localtop_t      *top,
252                                 const matrix               box,
253                                 const rvec                 x[],
254                                 rvec                       f[],
255                                 gmx::ForceWithVirial      *forceWithVirial,
256                                 tensor                     vir_force,
257                                 const t_mdatoms           *mdatoms,
258                                 const t_graph             *graph,
259                                 const t_forcerec          *fr,
260                                 const gmx_vsite_t         *vsite,
261                                 int                        flags)
262 {
263     if (fr->haveDirectVirialContributions)
264     {
265         rvec *fDirectVir = as_rvec_array(forceWithVirial->force_.data());
266
267         if (vsite)
268         {
269             /* Spread the mesh force on virtual sites to the other particles...
270              * This is parallellized. MPI communication is performed
271              * if the constructing atoms aren't local.
272              */
273             matrix virial = { { 0 } };
274             spread_vsite_f(vsite, x, fDirectVir, nullptr,
275                            (flags & GMX_FORCE_VIRIAL) != 0, virial,
276                            nrnb,
277                            &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
278             forceWithVirial->addVirialContribution(virial);
279         }
280
281         if (flags & GMX_FORCE_VIRIAL)
282         {
283             /* Now add the forces, this is local */
284             sum_forces(f, forceWithVirial->force_);
285
286             /* Add the direct virial contributions */
287             GMX_ASSERT(forceWithVirial->computeVirial_, "forceWithVirial should request virial computation when we request the virial");
288             m_add(vir_force, forceWithVirial->getVirial(), vir_force);
289
290             if (debug)
291             {
292                 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
293             }
294         }
295     }
296
297     if (fr->print_force >= 0)
298     {
299         print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
300     }
301 }
302
303 static void do_nb_verlet(t_forcerec                       *fr,
304                          const interaction_const_t        *ic,
305                          gmx_enerdata_t                   *enerd,
306                          const int                         flags,
307                          const Nbnxm::InteractionLocality  ilocality,
308                          const int                         clearF,
309                          const int64_t                     step,
310                          t_nrnb                           *nrnb,
311                          gmx_wallcycle_t                   wcycle)
312 {
313     if (!(flags & GMX_FORCE_NONBONDED))
314     {
315         /* skip non-bonded calculation */
316         return;
317     }
318
319     nonbonded_verlet_t *nbv  = fr->nbv.get();
320
321     /* GPU kernel launch overhead is already timed separately */
322     if (fr->cutoff_scheme != ecutsVERLET)
323     {
324         gmx_incons("Invalid cut-off scheme passed!");
325     }
326
327     if (!nbv->useGpu())
328     {
329         /* When dynamic pair-list  pruning is requested, we need to prune
330          * at nstlistPrune steps.
331          */
332         if (nbv->isDynamicPruningStepCpu(step))
333         {
334             /* Prune the pair-list beyond fr->ic->rlistPrune using
335              * the current coordinates of the atoms.
336              */
337             wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
338             nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
339             wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
340         }
341
342         wallcycle_sub_start(wcycle, ewcsNONBONDED);
343     }
344
345     nbv->dispatchNonbondedKernel(ilocality, *ic, flags, clearF, fr, enerd, nrnb);
346
347     if (!nbv->useGpu())
348     {
349         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
350     }
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                   int                            forceFlags,
471                   const gmx_edsam               *ed)
472 {
473     const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
474
475     return
476         ((computeForces && forceProviders->hasForceProvider()) ||         // forceProviders
477          (inputrec->bPull && pull_have_potential(inputrec->pull_work)) || // pull
478          inputrec->bRot ||                                                // enforced rotation
479          (ed != nullptr) ||                                               // flooding
480          (inputrec->bIMD && computeForces));                              // IMD
481 }
482
483 /*! \brief Compute forces and/or energies for special algorithms
484  *
485  * The intention is to collect all calls to algorithms that compute
486  * forces on local atoms only and that do not contribute to the local
487  * virial sum (but add their virial contribution separately).
488  * Eventually these should likely all become ForceProviders.
489  * Within this function the intention is to have algorithms that do
490  * global communication at the end, so global barriers within the MD loop
491  * are as close together as possible.
492  *
493  * \param[in]     fplog            The log file
494  * \param[in]     cr               The communication record
495  * \param[in]     inputrec         The input record
496  * \param[in]     awh              The Awh module (nullptr if none in use).
497  * \param[in]     enforcedRotation Enforced rotation module.
498  * \param[in]     step             The current MD step
499  * \param[in]     t                The current time
500  * \param[in,out] wcycle           Wallcycle accounting struct
501  * \param[in,out] forceProviders   Pointer to a list of force providers
502  * \param[in]     box              The unit cell
503  * \param[in]     x                The coordinates
504  * \param[in]     mdatoms          Per atom properties
505  * \param[in]     lambda           Array of free-energy lambda values
506  * \param[in]     forceFlags       Flags that tell whether we should compute forces/energies/virial
507  * \param[in,out] forceWithVirial  Force and virial buffers
508  * \param[in,out] enerd            Energy buffer
509  * \param[in,out] ed               Essential dynamics pointer
510  * \param[in]     bNS              Tells if we did neighbor searching this step, used for ED sampling
511  *
512  * \todo Remove bNS, which is used incorrectly.
513  * \todo Convert all other algorithms called here to ForceProviders.
514  */
515 static void
516 computeSpecialForces(FILE                          *fplog,
517                      const t_commrec               *cr,
518                      const t_inputrec              *inputrec,
519                      gmx::Awh                      *awh,
520                      gmx_enfrot                    *enforcedRotation,
521                      int64_t                        step,
522                      double                         t,
523                      gmx_wallcycle_t                wcycle,
524                      ForceProviders                *forceProviders,
525                      matrix                         box,
526                      gmx::ArrayRef<const gmx::RVec> x,
527                      const t_mdatoms               *mdatoms,
528                      real                          *lambda,
529                      int                            forceFlags,
530                      gmx::ForceWithVirial          *forceWithVirial,
531                      gmx_enerdata_t                *enerd,
532                      gmx_edsam                     *ed,
533                      gmx_bool                       bNS)
534 {
535     const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
536
537     /* NOTE: Currently all ForceProviders only provide forces.
538      *       When they also provide energies, remove this conditional.
539      */
540     if (computeForces)
541     {
542         gmx::ForceProviderInput  forceProviderInput(x, *mdatoms, t, box, *cr);
543         gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
544
545         /* Collect forces from modules */
546         forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
547     }
548
549     if (inputrec->bPull && pull_have_potential(inputrec->pull_work))
550     {
551         pull_potential_wrapper(cr, inputrec, box, x,
552                                forceWithVirial,
553                                mdatoms, enerd, lambda, t,
554                                wcycle);
555
556         if (awh)
557         {
558             enerd->term[F_COM_PULL] +=
559                 awh->applyBiasForcesAndUpdateBias(inputrec->ePBC, *mdatoms, box,
560                                                   forceWithVirial,
561                                                   t, step, wcycle, fplog);
562         }
563     }
564
565     rvec *f = as_rvec_array(forceWithVirial->force_.data());
566
567     /* Add the forces from enforced rotation potentials (if any) */
568     if (inputrec->bRot)
569     {
570         wallcycle_start(wcycle, ewcROTadd);
571         enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
572         wallcycle_stop(wcycle, ewcROTadd);
573     }
574
575     if (ed)
576     {
577         /* Note that since init_edsam() is called after the initialization
578          * of forcerec, edsam doesn't request the noVirSum force buffer.
579          * Thus if no other algorithm (e.g. PME) requires it, the forces
580          * here will contribute to the virial.
581          */
582         do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, bNS);
583     }
584
585     /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
586     if (inputrec->bIMD && computeForces)
587     {
588         IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
589     }
590 }
591
592 /*! \brief Launch the prepare_step and spread stages of PME GPU.
593  *
594  * \param[in]  pmedata       The PME structure
595  * \param[in]  box           The box matrix
596  * \param[in]  x             Coordinate array
597  * \param[in]  flags         Force flags
598  * \param[in]  pmeFlags      PME flags
599  * \param[in]  wcycle        The wallcycle structure
600  */
601 static inline void launchPmeGpuSpread(gmx_pme_t      *pmedata,
602                                       matrix          box,
603                                       rvec            x[],
604                                       int             flags,
605                                       int             pmeFlags,
606                                       gmx_wallcycle_t wcycle)
607 {
608     pme_gpu_prepare_computation(pmedata, (flags & GMX_FORCE_DYNAMICBOX) != 0, box, wcycle, pmeFlags);
609     pme_gpu_launch_spread(pmedata, x, wcycle);
610 }
611
612 /*! \brief Launch the FFT and gather stages of PME GPU
613  *
614  * This function only implements setting the output forces (no accumulation).
615  *
616  * \param[in]  pmedata        The PME structure
617  * \param[in]  wcycle         The wallcycle structure
618  */
619 static void launchPmeGpuFftAndGather(gmx_pme_t        *pmedata,
620                                      gmx_wallcycle_t   wcycle)
621 {
622     pme_gpu_launch_complex_transforms(pmedata, wcycle);
623     pme_gpu_launch_gather(pmedata, wcycle, PmeForceOutputHandling::Set);
624 }
625
626 /*! \brief
627  *  Polling wait for either of the PME or nonbonded GPU tasks.
628  *
629  * Instead of a static order in waiting for GPU tasks, this function
630  * polls checking which of the two tasks completes first, and does the
631  * associated force buffer reduction overlapped with the other task.
632  * By doing that, unlike static scheduling order, it can always overlap
633  * one of the reductions, regardless of the GPU task completion order.
634  *
635  * \param[in]     nbv              Nonbonded verlet structure
636  * \param[in,out] pmedata          PME module data
637  * \param[in,out] force            Force array to reduce task outputs into.
638  * \param[in,out] forceWithVirial  Force and virial buffers
639  * \param[in,out] fshift           Shift force output vector results are reduced into
640  * \param[in,out] enerd            Energy data structure results are reduced into
641  * \param[in]     flags            Force flags
642  * \param[in]     pmeFlags         PME flags
643  * \param[in]     haveOtherWork    Tells whether there is other work than non-bonded in the stream(s)
644  * \param[in]     wcycle           The wallcycle structure
645  */
646 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t                  *nbv,
647                                         gmx_pme_t                           *pmedata,
648                                         gmx::ArrayRefWithPadding<gmx::RVec> *force,
649                                         gmx::ForceWithVirial                *forceWithVirial,
650                                         rvec                                 fshift[],
651                                         gmx_enerdata_t                      *enerd,
652                                         int                                  flags,
653                                         int                                  pmeFlags,
654                                         bool                                 haveOtherWork,
655                                         gmx_wallcycle_t                      wcycle)
656 {
657     bool isPmeGpuDone = false;
658     bool isNbGpuDone  = false;
659
660
661     gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
662
663     while (!isPmeGpuDone || !isNbGpuDone)
664     {
665         if (!isPmeGpuDone)
666         {
667             GpuTaskCompletion completionType = (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
668             isPmeGpuDone = pme_gpu_try_finish_task(pmedata, pmeFlags, wcycle, forceWithVirial, enerd, completionType);
669         }
670
671         if (!isNbGpuDone)
672         {
673             GpuTaskCompletion completionType = (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
674             wallcycle_start_nocount(wcycle, ewcWAIT_GPU_NB_L);
675             isNbGpuDone = Nbnxm::gpu_try_finish_task(nbv->gpu_nbv,
676                                                      flags,
677                                                      Nbnxm::AtomLocality::Local,
678                                                      haveOtherWork,
679                                                      enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
680                                                      fshift, completionType);
681             wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
682             // To get the call count right, when the task finished we
683             // issue a start/stop.
684             // TODO: move the ewcWAIT_GPU_NB_L cycle counting into nbnxn_gpu_try_finish_task()
685             // and ewcNB_XF_BUF_OPS counting into nbnxn_atomdata_add_nbat_f_to_f().
686             if (isNbGpuDone)
687             {
688                 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
689                 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
690
691                 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
692                                               as_rvec_array(force->unpaddedArrayRef().data()), wcycle);
693             }
694         }
695     }
696 }
697
698 /*! \brief Hack structure with force ouput buffers for do_force for the home atoms for this domain */
699 struct ForceOutputs
700 {
701     //! Constructor
702     ForceOutputs(rvec *f, gmx::ForceWithVirial const forceWithVirial) :
703         f(f),
704         forceWithVirial(forceWithVirial) {}
705
706     //! Force output buffer used by legacy modules (without SIMD padding)
707     rvec                 *const f;
708     //! Force with direct virial contribution (if there are any; without SIMD padding)
709     gmx::ForceWithVirial        forceWithVirial;
710 };
711
712 /*! \brief Set up the different force buffers; also does clearing.
713  *
714  * \param[in] fr        force record pointer
715  * \param[in] inputrec  input record
716  * \param[in] force     force array
717  * \param[in] bDoForces True if force are computed this step
718  * \param[in] doVirial  True if virial is computed this step
719  * \param[out] wcycle   wallcycle recording structure
720  *
721  * \returns             Cleared force output structure
722  */
723 static ForceOutputs
724 setupForceOutputs(const t_forcerec                    *fr,
725                   const t_inputrec                    &inputrec,
726                   gmx::ArrayRefWithPadding<gmx::RVec>  force,
727                   const bool                           bDoForces,
728                   const bool                           doVirial,
729                   gmx_wallcycle_t                      wcycle)
730 {
731     wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
732
733     /* Temporary solution until all routines take PaddedRVecVector */
734     rvec *const f = as_rvec_array(force.unpaddedArrayRef().data());
735     if (bDoForces)
736     {
737         /* Clear the short- and long-range forces */
738         clear_rvecs_omp(fr->natoms_force_constr, f);
739     }
740
741     /* If we need to compute the virial, we might need a separate
742      * force buffer for algorithms for which the virial is calculated
743      * directly, such as PME. Otherwise, forceWithVirial uses the
744      * the same force (f in legacy calls) buffer as other algorithms.
745      */
746     const bool useSeparateForceWithVirialBuffer = (bDoForces && (doVirial && fr->haveDirectVirialContributions));
747
748
749     /* forceWithVirial uses the local atom range only */
750     gmx::ForceWithVirial forceWithVirial (useSeparateForceWithVirialBuffer ?
751                                           *fr->forceBufferForDirectVirialContributions : force.unpaddedArrayRef(),
752                                           doVirial);
753
754     if (useSeparateForceWithVirialBuffer)
755     {
756         /* TODO: update comment
757          * We only compute forces on local atoms. Note that vsites can
758          * spread to non-local atoms, but that part of the buffer is
759          * cleared separately in the vsite spreading code.
760          */
761         clear_rvecs_omp(forceWithVirial.force_.size(), as_rvec_array(forceWithVirial.force_.data()));
762     }
763
764     if (inputrec.bPull && pull_have_constraint(inputrec.pull_work))
765     {
766         clear_pull_forces(inputrec.pull_work);
767     }
768
769     wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
770
771     return ForceOutputs(f, forceWithVirial);
772 }
773
774
775 /*! \brief Set up flags that indicate what type of work is there to compute.
776  *
777  * Currently we only update it at search steps,
778  * but some properties may change more frequently (e.g. virial/non-virial step),
779  * so when including those either the frequency of update (per-step) or the scope
780  * of a flag will change (i.e. a set of flags for nstlist steps).
781  *
782  */
783 static void
784 setupForceWorkload(gmx::PpForceWorkload *forceWork,
785                    const t_inputrec     *inputrec,
786                    const t_forcerec     *fr,
787                    const gmx_edsam      *ed,
788                    const t_idef         &idef,
789                    const t_fcdata       *fcd,
790                    const int             forceFlags
791                    )
792 {
793     forceWork->haveSpecialForces      = haveSpecialForces(inputrec, fr->forceProviders, forceFlags, ed);
794     forceWork->haveCpuBondedWork      = haveCpuBondeds(*fr);
795     forceWork->haveGpuBondedWork      = ((fr->gpuBonded != nullptr) && fr->gpuBonded->haveInteractions());
796     forceWork->haveRestraintsWork     = havePositionRestraints(idef, *fcd);
797     forceWork->haveCpuListedForceWork = haveCpuListedForces(*fr, idef, *fcd);
798 }
799
800 static void do_force_cutsVERLET(FILE *fplog,
801                                 const t_commrec *cr,
802                                 const gmx_multisim_t *ms,
803                                 const t_inputrec *inputrec,
804                                 gmx::Awh *awh,
805                                 gmx_enfrot *enforcedRotation,
806                                 int64_t step,
807                                 t_nrnb *nrnb,
808                                 gmx_wallcycle_t wcycle,
809                                 const gmx_localtop_t *top,
810                                 const SimulationGroups * /* groups */,
811                                 matrix box, gmx::ArrayRefWithPadding<gmx::RVec> x,
812                                 history_t *hist,
813                                 gmx::ArrayRefWithPadding<gmx::RVec> force,
814                                 tensor vir_force,
815                                 const t_mdatoms *mdatoms,
816                                 gmx_enerdata_t *enerd, t_fcdata *fcd,
817                                 real *lambda,
818                                 t_graph *graph,
819                                 t_forcerec *fr,
820                                 gmx::PpForceWorkload *ppForceWorkload,
821                                 interaction_const_t *ic,
822                                 const gmx_vsite_t *vsite,
823                                 rvec mu_tot,
824                                 double t,
825                                 gmx_edsam *ed,
826                                 const int flags,
827                                 const DDBalanceRegionHandler &ddBalanceRegionHandler)
828 {
829     int                 i, j;
830     double              mu[2*DIM];
831     gmx_bool            bStateChanged, bNS, bFillGrid, bCalcCGCM;
832     gmx_bool            bDoForces, bUseGPU, bUseOrEmulGPU;
833     nonbonded_verlet_t *nbv = fr->nbv.get();
834
835     bStateChanged = ((flags & GMX_FORCE_STATECHANGED) != 0);
836     bNS           = ((flags & GMX_FORCE_NS) != 0);
837     bFillGrid     = (bNS && bStateChanged);
838     bCalcCGCM     = (bFillGrid && !DOMAINDECOMP(cr));
839     bDoForces     = ((flags & GMX_FORCE_FORCES) != 0);
840     bUseGPU       = fr->nbv->useGpu();
841     bUseOrEmulGPU = bUseGPU || fr->nbv->emulateGpu();
842
843     const auto pmeRunMode = fr->pmedata ? pme_run_mode(fr->pmedata) : PmeRunMode::CPU;
844     // TODO slim this conditional down - inputrec and duty checks should mean the same in proper code!
845     const bool useGpuPme  = EEL_PME(fr->ic->eeltype) && thisRankHasDuty(cr, DUTY_PME) &&
846         ((pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed));
847     const int  pmeFlags = GMX_PME_SPREAD | GMX_PME_SOLVE |
848         ((flags & GMX_FORCE_VIRIAL) ? GMX_PME_CALC_ENER_VIR : 0) |
849         ((flags & GMX_FORCE_ENERGY) ? GMX_PME_CALC_ENER_VIR : 0) |
850         ((flags & GMX_FORCE_FORCES) ? GMX_PME_CALC_F : 0);
851
852     /* At a search step we need to start the first balancing region
853      * somewhere early inside the step after communication during domain
854      * decomposition (and not during the previous step as usual).
855      */
856     if (bNS)
857     {
858         ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
859     }
860
861     const int start  = 0;
862     const int homenr = mdatoms->homenr;
863
864     clear_mat(vir_force);
865
866     if (bStateChanged)
867     {
868         update_forcerec(fr, box);
869
870         if (inputrecNeedMutot(inputrec))
871         {
872             /* Calculate total (local) dipole moment in a temporary common array.
873              * This makes it possible to sum them over nodes faster.
874              */
875             calc_mu(start, homenr,
876                     x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
877                     mu, mu+DIM);
878         }
879     }
880
881     if (fr->ePBC != epbcNONE)
882     {
883         /* Compute shift vectors every step,
884          * because of pressure coupling or box deformation!
885          */
886         if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
887         {
888             calc_shifts(box, fr->shift_vec);
889         }
890
891         if (bCalcCGCM)
892         {
893             put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr), gmx_omp_nthreads_get(emntDefault));
894             inc_nrnb(nrnb, eNR_SHIFTX, homenr);
895         }
896         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
897         {
898             unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
899         }
900     }
901
902     nbnxn_atomdata_copy_shiftvec((flags & GMX_FORCE_DYNAMICBOX) != 0,
903                                  fr->shift_vec, nbv->nbat.get());
904
905 #if GMX_MPI
906     if (!thisRankHasDuty(cr, DUTY_PME))
907     {
908         /* Send particle coordinates to the pme nodes.
909          * Since this is only implemented for domain decomposition
910          * and domain decomposition does not use the graph,
911          * we do not need to worry about shifting.
912          */
913         gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
914                                  lambda[efptCOUL], lambda[efptVDW],
915                                  (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)) != 0,
916                                  step, wcycle);
917     }
918 #endif /* GMX_MPI */
919
920     if (useGpuPme)
921     {
922         launchPmeGpuSpread(fr->pmedata, box, as_rvec_array(x.unpaddedArrayRef().data()), flags, pmeFlags, wcycle);
923     }
924
925     /* do gridding for pair search */
926     if (bNS)
927     {
928         if (graph && bStateChanged)
929         {
930             /* Calculate intramolecular shift vectors to make molecules whole */
931             mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
932         }
933
934         // TODO
935         // - vzero is constant, do we need to pass it?
936         // - box_diag should be passed directly to nbnxn_put_on_grid
937         //
938         rvec vzero;
939         clear_rvec(vzero);
940
941         rvec box_diag;
942         box_diag[XX] = box[XX][XX];
943         box_diag[YY] = box[YY][YY];
944         box_diag[ZZ] = box[ZZ][ZZ];
945
946         wallcycle_start(wcycle, ewcNS);
947         if (!DOMAINDECOMP(cr))
948         {
949             wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
950             nbnxn_put_on_grid(nbv, box,
951                               0, vzero, box_diag,
952                               nullptr, 0, mdatoms->homenr, -1,
953                               fr->cginfo, x.unpaddedArrayRef(),
954                               0, nullptr);
955             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
956         }
957         else
958         {
959             wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
960             nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd),
961                                        fr->cginfo, x.unpaddedArrayRef());
962             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
963         }
964
965         nbv->setAtomProperties(*mdatoms, *fr->cginfo);
966
967         wallcycle_stop(wcycle, ewcNS);
968
969         /* initialize the GPU nbnxm atom data and bonded data structures */
970         if (bUseGPU)
971         {
972             wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
973
974             wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
975             Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
976             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
977
978             if (fr->gpuBonded)
979             {
980                 /* Now we put all atoms on the grid, we can assign bonded
981                  * interactions to the GPU, where the grid order is
982                  * needed. Also the xq, f and fshift device buffers have
983                  * been reallocated if needed, so the bonded code can
984                  * learn about them. */
985                 // TODO the xq, f, and fshift buffers are now shared
986                 // resources, so they should be maintained by a
987                 // higher-level object than the nb module.
988                 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
989                                                                       top->idef,
990                                                                       Nbnxm::gpu_get_xq(nbv->gpu_nbv),
991                                                                       Nbnxm::gpu_get_f(nbv->gpu_nbv),
992                                                                       Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
993             }
994             wallcycle_stop(wcycle, ewcLAUNCH_GPU);
995         }
996
997         // Need to run after the GPU-offload bonded interaction lists
998         // are set up to be able to determine whether there is bonded work.
999         setupForceWorkload(ppForceWorkload,
1000                            inputrec,
1001                            fr,
1002                            ed,
1003                            top->idef,
1004                            fcd,
1005                            flags);
1006     }
1007
1008     /* do local pair search */
1009     if (bNS)
1010     {
1011         // TODO: fuse this branch with the above bNS block
1012         wallcycle_start_nocount(wcycle, ewcNS);
1013         wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1014         /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1015         nbv->constructPairlist(Nbnxm::InteractionLocality::Local,
1016                                &top->excls, step, nrnb);
1017         wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1018         wallcycle_stop(wcycle, ewcNS);
1019     }
1020     else
1021     {
1022         nbv->setCoordinates(Nbnxm::AtomLocality::Local, false,
1023                             x.unpaddedArrayRef(), wcycle);
1024     }
1025
1026     if (bUseGPU)
1027     {
1028         ddBalanceRegionHandler.openBeforeForceComputationGpu();
1029
1030         wallcycle_start(wcycle, ewcLAUNCH_GPU);
1031
1032         wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1033         Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1034         Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1035                                   Nbnxm::AtomLocality::Local,
1036                                   ppForceWorkload->haveGpuBondedWork);
1037         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1038
1039         // bonded work not split into separate local and non-local, so with DD
1040         // we can only launch the kernel after non-local coordinates have been received.
1041         if (ppForceWorkload->haveGpuBondedWork && !havePPDomainDecomposition(cr))
1042         {
1043             wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1044             fr->gpuBonded->launchKernels(fr, flags, box);
1045             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1046         }
1047
1048         /* launch local nonbonded work on GPU */
1049         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1050         do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFNo,
1051                      step, nrnb, wcycle);
1052         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1053         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1054     }
1055
1056     if (useGpuPme)
1057     {
1058         // In PME GPU and mixed mode we launch FFT / gather after the
1059         // X copy/transform to allow overlap as well as after the GPU NB
1060         // launch to avoid FFT launch overhead hijacking the CPU and delaying
1061         // the nonbonded kernel.
1062         launchPmeGpuFftAndGather(fr->pmedata, wcycle);
1063     }
1064
1065     /* Communicate coordinates and sum dipole if necessary +
1066        do non-local pair search */
1067     if (havePPDomainDecomposition(cr))
1068     {
1069         if (bNS)
1070         {
1071             // TODO: fuse this branch with the above large bNS block
1072             wallcycle_start_nocount(wcycle, ewcNS);
1073             wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1074             /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1075             nbv->constructPairlist(Nbnxm::InteractionLocality::NonLocal,
1076                                    &top->excls, step, nrnb);
1077             wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1078             wallcycle_stop(wcycle, ewcNS);
1079         }
1080         else
1081         {
1082             dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1083
1084             nbv->setCoordinates(Nbnxm::AtomLocality::NonLocal, false,
1085                                 x.unpaddedArrayRef(), wcycle);
1086         }
1087
1088         if (bUseGPU)
1089         {
1090             wallcycle_start(wcycle, ewcLAUNCH_GPU);
1091
1092             /* launch non-local nonbonded tasks on GPU */
1093             wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1094             Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1095                                       Nbnxm::AtomLocality::NonLocal,
1096                                       ppForceWorkload->haveGpuBondedWork);
1097             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1098
1099             if (ppForceWorkload->haveGpuBondedWork)
1100             {
1101                 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1102                 fr->gpuBonded->launchKernels(fr, flags, box);
1103                 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1104             }
1105
1106             wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1107             do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1108                          step, nrnb, wcycle);
1109             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1110
1111             wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1112         }
1113     }
1114
1115     if (bUseGPU)
1116     {
1117         /* launch D2H copy-back F */
1118         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1119         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1120         if (havePPDomainDecomposition(cr))
1121         {
1122             Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1123                                       flags, Nbnxm::AtomLocality::NonLocal, ppForceWorkload->haveGpuBondedWork);
1124         }
1125         Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1126                                   flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork);
1127         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1128
1129         if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1130         {
1131             wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1132             fr->gpuBonded->launchEnergyTransfer();
1133             wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1134         }
1135         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1136     }
1137
1138     if (bStateChanged && inputrecNeedMutot(inputrec))
1139     {
1140         if (PAR(cr))
1141         {
1142             gmx_sumd(2*DIM, mu, cr);
1143
1144             ddBalanceRegionHandler.reopenRegionCpu();
1145         }
1146
1147         for (i = 0; i < 2; i++)
1148         {
1149             for (j = 0; j < DIM; j++)
1150             {
1151                 fr->mu_tot[i][j] = mu[i*DIM + j];
1152             }
1153         }
1154     }
1155     if (fr->efep == efepNO)
1156     {
1157         copy_rvec(fr->mu_tot[0], mu_tot);
1158     }
1159     else
1160     {
1161         for (j = 0; j < DIM; j++)
1162         {
1163             mu_tot[j] =
1164                 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1165                 lambda[efptCOUL]*fr->mu_tot[1][j];
1166         }
1167     }
1168
1169     /* Reset energies */
1170     reset_enerdata(enerd);
1171     clear_rvecs(SHIFTS, fr->fshift);
1172
1173     if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1174     {
1175         wallcycle_start(wcycle, ewcPPDURINGPME);
1176         dd_force_flop_start(cr->dd, nrnb);
1177     }
1178
1179     if (inputrec->bRot)
1180     {
1181         wallcycle_start(wcycle, ewcROT);
1182         do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, bNS);
1183         wallcycle_stop(wcycle, ewcROT);
1184     }
1185
1186     /* Start the force cycle counter.
1187      * Note that a different counter is used for dynamic load balancing.
1188      */
1189     wallcycle_start(wcycle, ewcFORCE);
1190
1191     // set up and clear force outputs
1192     struct ForceOutputs forceOut = setupForceOutputs(fr, *inputrec, force, bDoForces, ((flags & GMX_FORCE_VIRIAL) != 0), wcycle);
1193
1194     /* We calculate the non-bonded forces, when done on the CPU, here.
1195      * We do this before calling do_force_lowlevel, because in that
1196      * function, the listed forces are calculated before PME, which
1197      * does communication.  With this order, non-bonded and listed
1198      * force calculation imbalance can be balanced out by the domain
1199      * decomposition load balancing.
1200      */
1201
1202     if (!bUseOrEmulGPU)
1203     {
1204         do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFYes,
1205                      step, nrnb, wcycle);
1206     }
1207
1208     if (fr->efep != efepNO)
1209     {
1210         /* Calculate the local and non-local free energy interactions here.
1211          * Happens here on the CPU both with and without GPU.
1212          */
1213         wallcycle_sub_start(wcycle, ewcsNONBONDED);
1214         nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::Local,
1215                                       fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1216                                       inputrec->fepvals, lambda,
1217                                       enerd, flags, nrnb);
1218
1219         if (havePPDomainDecomposition(cr))
1220         {
1221             nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::NonLocal,
1222                                           fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1223                                           inputrec->fepvals, lambda,
1224                                           enerd, flags, nrnb);
1225         }
1226         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
1227     }
1228
1229     if (!bUseOrEmulGPU)
1230     {
1231         if (havePPDomainDecomposition(cr))
1232         {
1233             do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1234                          step, nrnb, wcycle);
1235         }
1236
1237         /* Add all the non-bonded force to the normal force array.
1238          * This can be split into a local and a non-local part when overlapping
1239          * communication with calculation with domain decomposition.
1240          */
1241         wallcycle_stop(wcycle, ewcFORCE);
1242
1243         nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::All, forceOut.f, wcycle);
1244
1245         wallcycle_start_nocount(wcycle, ewcFORCE);
1246
1247         /* If there are multiple fshift output buffers we need to reduce them */
1248         if (flags & GMX_FORCE_VIRIAL)
1249         {
1250             /* This is not in a subcounter because it takes a
1251                negligible and constant-sized amount of time */
1252             nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->nbat.get(),
1253                                                      fr->fshift);
1254         }
1255     }
1256
1257     /* update QMMMrec, if necessary */
1258     if (fr->bQMMM)
1259     {
1260         update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1261     }
1262
1263     /* Compute the bonded and non-bonded energies and optionally forces */
1264     do_force_lowlevel(fr, inputrec, &(top->idef),
1265                       cr, ms, nrnb, wcycle, mdatoms,
1266                       as_rvec_array(x.unpaddedArrayRef().data()), hist, forceOut.f, &forceOut.forceWithVirial, enerd, fcd,
1267                       box, inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1268                       flags,
1269                       ddBalanceRegionHandler);
1270
1271     wallcycle_stop(wcycle, ewcFORCE);
1272
1273     computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1274                          step, t, wcycle,
1275                          fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda,
1276                          flags, &forceOut.forceWithVirial, enerd,
1277                          ed, bNS);
1278
1279     // Will store the amount of cycles spent waiting for the GPU that
1280     // will be later used in the DLB accounting.
1281     float cycles_wait_gpu = 0;
1282     if (bUseOrEmulGPU)
1283     {
1284         /* wait for non-local forces (or calculate in emulation mode) */
1285         if (havePPDomainDecomposition(cr))
1286         {
1287             if (bUseGPU)
1288             {
1289                 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1290                 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1291                                             flags, Nbnxm::AtomLocality::NonLocal,
1292                                             ppForceWorkload->haveGpuBondedWork,
1293                                             enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1294                                             fr->fshift);
1295                 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1296             }
1297             else
1298             {
1299                 wallcycle_start_nocount(wcycle, ewcFORCE);
1300                 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFYes,
1301                              step, nrnb, wcycle);
1302                 wallcycle_stop(wcycle, ewcFORCE);
1303             }
1304
1305             nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::NonLocal,
1306                                           forceOut.f, wcycle);
1307         }
1308     }
1309
1310     if (havePPDomainDecomposition(cr))
1311     {
1312         /* We are done with the CPU compute.
1313          * We will now communicate the non-local forces.
1314          * If we use a GPU this will overlap with GPU work, so in that case
1315          * we do not close the DD force balancing region here.
1316          */
1317         ddBalanceRegionHandler.closeAfterForceComputationCpu();
1318
1319         if (bDoForces)
1320         {
1321             dd_move_f(cr->dd, force.unpaddedArrayRef(), fr->fshift, wcycle);
1322         }
1323     }
1324
1325     // With both nonbonded and PME offloaded a GPU on the same rank, we use
1326     // an alternating wait/reduction scheme.
1327     bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPme && bUseGPU && !DOMAINDECOMP(cr));
1328     if (alternateGpuWait)
1329     {
1330         alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &force, &forceOut.forceWithVirial, fr->fshift, enerd,
1331                                     flags, pmeFlags, ppForceWorkload->haveGpuBondedWork, wcycle);
1332     }
1333
1334     if (!alternateGpuWait && useGpuPme)
1335     {
1336         pme_gpu_wait_and_reduce(fr->pmedata, pmeFlags, wcycle, &forceOut.forceWithVirial, enerd);
1337     }
1338
1339     /* Wait for local GPU NB outputs on the non-alternating wait path */
1340     if (!alternateGpuWait && bUseGPU)
1341     {
1342         /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1343          * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1344          * but even with a step of 0.1 ms the difference is less than 1%
1345          * of the step time.
1346          */
1347         const float gpuWaitApiOverheadMargin = 2e6f; /* cycles */
1348
1349         wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1350         Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1351                                     flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork,
1352                                     enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1353                                     fr->fshift);
1354         float cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1355
1356         if (ddBalanceRegionHandler.useBalancingRegion())
1357         {
1358             DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1359             if (bDoForces && cycles_tmp <= gpuWaitApiOverheadMargin)
1360             {
1361                 /* We measured few cycles, it could be that the kernel
1362                  * and transfer finished earlier and there was no actual
1363                  * wait time, only API call overhead.
1364                  * Then the actual time could be anywhere between 0 and
1365                  * cycles_wait_est. We will use half of cycles_wait_est.
1366                  */
1367                 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1368             }
1369             ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1370         }
1371     }
1372
1373     if (fr->nbv->emulateGpu())
1374     {
1375         // NOTE: emulation kernel is not included in the balancing region,
1376         // but emulation mode does not target performance anyway
1377         wallcycle_start_nocount(wcycle, ewcFORCE);
1378         do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local,
1379                      DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1380                      step, nrnb, wcycle);
1381         wallcycle_stop(wcycle, ewcFORCE);
1382     }
1383
1384     if (useGpuPme)
1385     {
1386         pme_gpu_reinit_computation(fr->pmedata, wcycle);
1387     }
1388
1389     if (bUseGPU)
1390     {
1391         /* now clear the GPU outputs while we finish the step on the CPU */
1392         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1393         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1394         Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, flags);
1395
1396         if (nbv->isDynamicPruningStepGpu(step))
1397         {
1398             nbv->dispatchPruneKernelGpu(step);
1399         }
1400         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1401         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1402     }
1403
1404     if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1405     {
1406         wallcycle_start(wcycle, ewcWAIT_GPU_BONDED);
1407         // in principle this should be included in the DD balancing region,
1408         // but generally it is infrequent so we'll omit it for the sake of
1409         // simpler code
1410         fr->gpuBonded->accumulateEnergyTerms(enerd);
1411         wallcycle_stop(wcycle, ewcWAIT_GPU_BONDED);
1412
1413         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1414         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1415         fr->gpuBonded->clearEnergies();
1416         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1417         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1418     }
1419
1420     /* Do the nonbonded GPU (or emulation) force buffer reduction
1421      * on the non-alternating path. */
1422     if (bUseOrEmulGPU && !alternateGpuWait)
1423     {
1424         nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
1425                                       forceOut.f, wcycle);
1426     }
1427     if (DOMAINDECOMP(cr))
1428     {
1429         dd_force_flop_stop(cr->dd, nrnb);
1430     }
1431
1432     if (bDoForces)
1433     {
1434         /* If we have NoVirSum forces, but we do not calculate the virial,
1435          * we sum fr->f_novirsum=forceOut.f later.
1436          */
1437         if (vsite && !(fr->haveDirectVirialContributions && !(flags & GMX_FORCE_VIRIAL)))
1438         {
1439             spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, fr->fshift, FALSE, nullptr, nrnb,
1440                            &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1441         }
1442
1443         if (flags & GMX_FORCE_VIRIAL)
1444         {
1445             /* Calculation of the virial must be done after vsites! */
1446             calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f,
1447                         vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1448         }
1449     }
1450
1451     if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1452     {
1453         /* In case of node-splitting, the PP nodes receive the long-range
1454          * forces, virial and energy from the PME nodes here.
1455          */
1456         pme_receive_force_ener(cr, &forceOut.forceWithVirial, enerd, wcycle);
1457     }
1458
1459     if (bDoForces)
1460     {
1461         post_process_forces(cr, step, nrnb, wcycle,
1462                             top, box, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, &forceOut.forceWithVirial,
1463                             vir_force, mdatoms, graph, fr, vsite,
1464                             flags);
1465     }
1466
1467     if (flags & GMX_FORCE_ENERGY)
1468     {
1469         /* Sum the potential energy terms from group contributions */
1470         sum_epot(&(enerd->grpp), enerd->term);
1471
1472         if (!EI_TPI(inputrec->eI))
1473         {
1474             checkPotentialEnergyValidity(step, *enerd, *inputrec);
1475         }
1476     }
1477 }
1478
1479 static void do_force_cutsGROUP(FILE *fplog,
1480                                const t_commrec *cr,
1481                                const gmx_multisim_t *ms,
1482                                const t_inputrec *inputrec,
1483                                gmx::Awh *awh,
1484                                gmx_enfrot *enforcedRotation,
1485                                int64_t step,
1486                                t_nrnb *nrnb,
1487                                gmx_wallcycle_t wcycle,
1488                                gmx_localtop_t *top,
1489                                const SimulationGroups *groups,
1490                                matrix box, gmx::ArrayRefWithPadding<gmx::RVec> x,
1491                                history_t *hist,
1492                                gmx::ArrayRefWithPadding<gmx::RVec> force,
1493                                tensor vir_force,
1494                                const t_mdatoms *mdatoms,
1495                                gmx_enerdata_t *enerd,
1496                                t_fcdata *fcd,
1497                                real *lambda,
1498                                t_graph *graph,
1499                                t_forcerec *fr,
1500                                const gmx_vsite_t *vsite,
1501                                rvec mu_tot,
1502                                double t,
1503                                gmx_edsam *ed,
1504                                int flags,
1505                                const DDBalanceRegionHandler &ddBalanceRegionHandler)
1506 {
1507     int        cg0, cg1, i, j;
1508     double     mu[2*DIM];
1509     gmx_bool   bStateChanged, bNS, bFillGrid, bCalcCGCM;
1510     gmx_bool   bDoForces;
1511
1512     const int  start  = 0;
1513     const int  homenr = mdatoms->homenr;
1514
1515     clear_mat(vir_force);
1516
1517     cg0 = 0;
1518     if (DOMAINDECOMP(cr))
1519     {
1520         cg1 = cr->dd->globalAtomGroupIndices.size();
1521     }
1522     else
1523     {
1524         cg1 = top->cgs.nr;
1525     }
1526     if (fr->n_tpi > 0)
1527     {
1528         cg1--;
1529     }
1530
1531     bStateChanged  = ((flags & GMX_FORCE_STATECHANGED) != 0);
1532     bNS            = ((flags & GMX_FORCE_NS) != 0);
1533     /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1534     bFillGrid      = (bNS && bStateChanged);
1535     bCalcCGCM      = (bFillGrid && !DOMAINDECOMP(cr));
1536     bDoForces      = ((flags & GMX_FORCE_FORCES) != 0);
1537
1538     if (bStateChanged)
1539     {
1540         update_forcerec(fr, box);
1541
1542         if (inputrecNeedMutot(inputrec))
1543         {
1544             /* Calculate total (local) dipole moment in a temporary common array.
1545              * This makes it possible to sum them over nodes faster.
1546              */
1547             calc_mu(start, homenr,
1548                     x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1549                     mu, mu+DIM);
1550         }
1551     }
1552
1553     if (fr->ePBC != epbcNONE)
1554     {
1555         /* Compute shift vectors every step,
1556          * because of pressure coupling or box deformation!
1557          */
1558         if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1559         {
1560             calc_shifts(box, fr->shift_vec);
1561         }
1562
1563         if (bCalcCGCM)
1564         {
1565             put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1566                                      &(top->cgs), as_rvec_array(x.unpaddedArrayRef().data()), fr->cg_cm);
1567             inc_nrnb(nrnb, eNR_CGCM, homenr);
1568             inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1569         }
1570         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1571         {
1572             unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
1573         }
1574     }
1575     else if (bCalcCGCM)
1576     {
1577         calc_cgcm(fplog, cg0, cg1, &(top->cgs), as_rvec_array(x.unpaddedArrayRef().data()), fr->cg_cm);
1578         inc_nrnb(nrnb, eNR_CGCM, homenr);
1579     }
1580
1581     if (bCalcCGCM && gmx_debug_at)
1582     {
1583         pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1584     }
1585
1586 #if GMX_MPI
1587     if (!thisRankHasDuty(cr, DUTY_PME))
1588     {
1589         /* Send particle coordinates to the pme nodes.
1590          * Since this is only implemented for domain decomposition
1591          * and domain decomposition does not use the graph,
1592          * we do not need to worry about shifting.
1593          */
1594         gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
1595                                  lambda[efptCOUL], lambda[efptVDW],
1596                                  (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)) != 0,
1597                                  step, wcycle);
1598     }
1599 #endif /* GMX_MPI */
1600
1601     /* Communicate coordinates and sum dipole if necessary */
1602     if (DOMAINDECOMP(cr))
1603     {
1604         dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1605
1606         /* No GPU support, no move_x overlap, so reopen the balance region here */
1607         ddBalanceRegionHandler.reopenRegionCpu();
1608     }
1609
1610     if (inputrecNeedMutot(inputrec))
1611     {
1612         if (bStateChanged)
1613         {
1614             if (PAR(cr))
1615             {
1616                 gmx_sumd(2*DIM, mu, cr);
1617
1618                 ddBalanceRegionHandler.reopenRegionCpu();
1619             }
1620             for (i = 0; i < 2; i++)
1621             {
1622                 for (j = 0; j < DIM; j++)
1623                 {
1624                     fr->mu_tot[i][j] = mu[i*DIM + j];
1625                 }
1626             }
1627         }
1628         if (fr->efep == efepNO)
1629         {
1630             copy_rvec(fr->mu_tot[0], mu_tot);
1631         }
1632         else
1633         {
1634             for (j = 0; j < DIM; j++)
1635             {
1636                 mu_tot[j] =
1637                     (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1638             }
1639         }
1640     }
1641
1642     /* Reset energies */
1643     reset_enerdata(enerd);
1644     clear_rvecs(SHIFTS, fr->fshift);
1645
1646     if (bNS)
1647     {
1648         wallcycle_start(wcycle, ewcNS);
1649
1650         if (graph && bStateChanged)
1651         {
1652             /* Calculate intramolecular shift vectors to make molecules whole */
1653             mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
1654         }
1655
1656         /* Do the actual neighbour searching */
1657         ns(fplog, fr, box,
1658            groups, top, mdatoms,
1659            cr, nrnb, bFillGrid);
1660
1661         wallcycle_stop(wcycle, ewcNS);
1662     }
1663
1664     if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1665     {
1666         wallcycle_start(wcycle, ewcPPDURINGPME);
1667         dd_force_flop_start(cr->dd, nrnb);
1668     }
1669
1670     if (inputrec->bRot)
1671     {
1672         wallcycle_start(wcycle, ewcROT);
1673         do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, bNS);
1674         wallcycle_stop(wcycle, ewcROT);
1675     }
1676
1677
1678     /* Start the force cycle counter.
1679      * Note that a different counter is used for dynamic load balancing.
1680      */
1681     wallcycle_start(wcycle, ewcFORCE);
1682
1683     // set up and clear force outputs
1684     struct ForceOutputs forceOut = setupForceOutputs(fr, *inputrec, force, bDoForces, ((flags & GMX_FORCE_VIRIAL) != 0), wcycle);
1685
1686     if (inputrec->bPull && pull_have_constraint(inputrec->pull_work))
1687     {
1688         clear_pull_forces(inputrec->pull_work);
1689     }
1690
1691     /* update QMMMrec, if necessary */
1692     if (fr->bQMMM)
1693     {
1694         update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1695     }
1696
1697     /* Compute the bonded and non-bonded energies and optionally forces */
1698     do_force_lowlevel(fr, inputrec, &(top->idef),
1699                       cr, ms, nrnb, wcycle, mdatoms,
1700                       as_rvec_array(x.unpaddedArrayRef().data()), hist, forceOut.f, &forceOut.forceWithVirial, enerd, fcd,
1701                       box, inputrec->fepvals, lambda,
1702                       graph, &(top->excls), fr->mu_tot,
1703                       flags,
1704                       ddBalanceRegionHandler);
1705
1706     wallcycle_stop(wcycle, ewcFORCE);
1707
1708     if (DOMAINDECOMP(cr))
1709     {
1710         dd_force_flop_stop(cr->dd, nrnb);
1711
1712         ddBalanceRegionHandler.closeAfterForceComputationCpu();
1713     }
1714
1715     computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1716                          step, t, wcycle,
1717                          fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda,
1718                          flags, &forceOut.forceWithVirial, enerd,
1719                          ed, bNS);
1720
1721     if (bDoForces)
1722     {
1723         /* Communicate the forces */
1724         if (DOMAINDECOMP(cr))
1725         {
1726             dd_move_f(cr->dd, force.unpaddedArrayRef(), fr->fshift, wcycle);
1727             /* Do we need to communicate the separate force array
1728              * for terms that do not contribute to the single sum virial?
1729              * Position restraints and electric fields do not introduce
1730              * inter-cg forces, only full electrostatics methods do.
1731              * When we do not calculate the virial, fr->f_novirsum = forceOut.f,
1732              * so we have already communicated these forces.
1733              */
1734             if (EEL_FULL(fr->ic->eeltype) && cr->dd->n_intercg_excl &&
1735                 (flags & GMX_FORCE_VIRIAL))
1736             {
1737                 dd_move_f(cr->dd, forceOut.forceWithVirial.force_, nullptr, wcycle);
1738             }
1739         }
1740
1741         /* If we have NoVirSum forces, but we do not calculate the virial,
1742          * we sum fr->f_novirsum=forceOut.f later.
1743          */
1744         if (vsite && !(fr->haveDirectVirialContributions && !(flags & GMX_FORCE_VIRIAL)))
1745         {
1746             spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, fr->fshift, FALSE, nullptr, nrnb,
1747                            &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1748         }
1749
1750         if (flags & GMX_FORCE_VIRIAL)
1751         {
1752             /* Calculation of the virial must be done after vsites! */
1753             calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f,
1754                         vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1755         }
1756     }
1757
1758     if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1759     {
1760         /* In case of node-splitting, the PP nodes receive the long-range
1761          * forces, virial and energy from the PME nodes here.
1762          */
1763         pme_receive_force_ener(cr, &forceOut.forceWithVirial, enerd, wcycle);
1764     }
1765
1766     if (bDoForces)
1767     {
1768         post_process_forces(cr, step, nrnb, wcycle,
1769                             top, box, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, &forceOut.forceWithVirial,
1770                             vir_force, mdatoms, graph, fr, vsite,
1771                             flags);
1772     }
1773
1774     if (flags & GMX_FORCE_ENERGY)
1775     {
1776         /* Sum the potential energy terms from group contributions */
1777         sum_epot(&(enerd->grpp), enerd->term);
1778
1779         if (!EI_TPI(inputrec->eI))
1780         {
1781             checkPotentialEnergyValidity(step, *enerd, *inputrec);
1782         }
1783     }
1784
1785 }
1786
1787 void do_force(FILE                                            *fplog,
1788               const t_commrec                                 *cr,
1789               const gmx_multisim_t                            *ms,
1790               const t_inputrec                                *inputrec,
1791               gmx::Awh                                        *awh,
1792               gmx_enfrot                                      *enforcedRotation,
1793               int64_t                                          step,
1794               t_nrnb                                          *nrnb,
1795               gmx_wallcycle_t                                  wcycle,
1796               gmx_localtop_t                                  *top,
1797               const SimulationGroups                          *groups,
1798               matrix                                           box,
1799               gmx::ArrayRefWithPadding<gmx::RVec>              x,     //NOLINT(performance-unnecessary-value-param)
1800               history_t                                       *hist,
1801               gmx::ArrayRefWithPadding<gmx::RVec>              force, //NOLINT(performance-unnecessary-value-param)
1802               tensor                                           vir_force,
1803               const t_mdatoms                                 *mdatoms,
1804               gmx_enerdata_t                                  *enerd,
1805               t_fcdata                                        *fcd,
1806               gmx::ArrayRef<real>                              lambda,
1807               t_graph                                         *graph,
1808               t_forcerec                                      *fr,
1809               gmx::PpForceWorkload                            *ppForceWorkload,
1810               const gmx_vsite_t                               *vsite,
1811               rvec                                             mu_tot,
1812               double                                           t,
1813               gmx_edsam                                       *ed,
1814               int                                              flags,
1815               const DDBalanceRegionHandler                    &ddBalanceRegionHandler)
1816 {
1817     /* modify force flag if not doing nonbonded */
1818     if (!fr->bNonbonded)
1819     {
1820         flags &= ~GMX_FORCE_NONBONDED;
1821     }
1822
1823     switch (inputrec->cutoff_scheme)
1824     {
1825         case ecutsVERLET:
1826             do_force_cutsVERLET(fplog, cr, ms, inputrec,
1827                                 awh, enforcedRotation, step, nrnb, wcycle,
1828                                 top,
1829                                 groups,
1830                                 box, x, hist,
1831                                 force, vir_force,
1832                                 mdatoms,
1833                                 enerd, fcd,
1834                                 lambda.data(), graph,
1835                                 fr,
1836                                 ppForceWorkload,
1837                                 fr->ic,
1838                                 vsite, mu_tot,
1839                                 t, ed,
1840                                 flags,
1841                                 ddBalanceRegionHandler);
1842             break;
1843         case ecutsGROUP:
1844             do_force_cutsGROUP(fplog, cr, ms, inputrec,
1845                                awh, enforcedRotation, step, nrnb, wcycle,
1846                                top,
1847                                groups,
1848                                box, x, hist,
1849                                force, vir_force,
1850                                mdatoms,
1851                                enerd, fcd,
1852                                lambda.data(), graph,
1853                                fr, vsite, mu_tot,
1854                                t, ed,
1855                                flags,
1856                                ddBalanceRegionHandler);
1857             break;
1858         default:
1859             gmx_incons("Invalid cut-off scheme passed!");
1860     }
1861
1862     /* In case we don't have constraints and are using GPUs, the next balancing
1863      * region starts here.
1864      * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1865      * virial calculation and COM pulling, is not thus not included in
1866      * the balance timing, which is ok as most tasks do communication.
1867      */
1868     ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);
1869 }