2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013-2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
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.
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.
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.
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.
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.
48 #include "gromacs/awh/awh.h"
49 #include "gromacs/domdec/dlbtiming.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/domdec/gpuhaloexchange.h"
53 #include "gromacs/domdec/partition.h"
54 #include "gromacs/essentialdynamics/edsam.h"
55 #include "gromacs/ewald/pme.h"
56 #include "gromacs/ewald/pme_pp.h"
57 #include "gromacs/ewald/pme_pp_comm_gpu.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
60 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
61 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
62 #include "gromacs/gpu_utils/gpu_utils.h"
63 #include "gromacs/imd/imd.h"
64 #include "gromacs/listed_forces/disre.h"
65 #include "gromacs/listed_forces/gpubonded.h"
66 #include "gromacs/listed_forces/listed_forces.h"
67 #include "gromacs/listed_forces/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/force_flags.h"
80 #include "gromacs/mdlib/forcerec.h"
81 #include "gromacs/mdlib/gmx_omp_nthreads.h"
82 #include "gromacs/mdlib/update.h"
83 #include "gromacs/mdlib/vsite.h"
84 #include "gromacs/mdlib/wholemoleculetransform.h"
85 #include "gromacs/mdtypes/commrec.h"
86 #include "gromacs/mdtypes/enerdata.h"
87 #include "gromacs/mdtypes/forceoutput.h"
88 #include "gromacs/mdtypes/forcerec.h"
89 #include "gromacs/mdtypes/iforceprovider.h"
90 #include "gromacs/mdtypes/inputrec.h"
91 #include "gromacs/mdtypes/md_enums.h"
92 #include "gromacs/mdtypes/mdatom.h"
93 #include "gromacs/mdtypes/simulation_workload.h"
94 #include "gromacs/mdtypes/state.h"
95 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
96 #include "gromacs/nbnxm/gpu_data_mgmt.h"
97 #include "gromacs/nbnxm/nbnxm.h"
98 #include "gromacs/nbnxm/nbnxm_gpu.h"
99 #include "gromacs/pbcutil/ishift.h"
100 #include "gromacs/pbcutil/pbc.h"
101 #include "gromacs/pulling/pull.h"
102 #include "gromacs/pulling/pull_rotation.h"
103 #include "gromacs/timing/cyclecounter.h"
104 #include "gromacs/timing/gpu_timing.h"
105 #include "gromacs/timing/wallcycle.h"
106 #include "gromacs/timing/wallcyclereporting.h"
107 #include "gromacs/timing/walltime_accounting.h"
108 #include "gromacs/topology/topology.h"
109 #include "gromacs/utility/arrayref.h"
110 #include "gromacs/utility/basedefinitions.h"
111 #include "gromacs/utility/cstringutil.h"
112 #include "gromacs/utility/exceptions.h"
113 #include "gromacs/utility/fatalerror.h"
114 #include "gromacs/utility/fixedcapacityvector.h"
115 #include "gromacs/utility/gmxassert.h"
116 #include "gromacs/utility/gmxmpi.h"
117 #include "gromacs/utility/logger.h"
118 #include "gromacs/utility/smalloc.h"
119 #include "gromacs/utility/strconvert.h"
120 #include "gromacs/utility/sysinfo.h"
123 using gmx::AtomLocality;
124 using gmx::DomainLifetimeWorkload;
125 using gmx::ForceOutputs;
126 using gmx::InteractionLocality;
128 using gmx::SimulationWorkload;
129 using gmx::StepWorkload;
131 // TODO: this environment variable allows us to verify before release
132 // that on less common architectures the total cost of polling is not larger than
133 // a blocking wait (so polling does not introduce overhead when the static
134 // PME-first ordering would suffice).
135 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
137 static void sum_forces(rvec f[], gmx::ArrayRef<const gmx::RVec> forceToAdd)
139 const int end = forceToAdd.size();
141 int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
142 #pragma omp parallel for num_threads(nt) schedule(static)
143 for (int i = 0; i < end; i++)
145 rvec_inc(f[i], forceToAdd[i]);
149 static void calc_virial(int start,
152 const gmx::ForceWithShiftForces& forceWithShiftForces,
156 const t_forcerec* fr,
159 /* The short-range virial from surrounding boxes */
160 const rvec* fshift = as_rvec_array(forceWithShiftForces.shiftForces().data());
161 calc_vir(SHIFTS, fr->shift_vec, fshift, vir_part, pbcType == PbcType::Screw, box);
162 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
164 /* Calculate partial virial, for local atoms only, based on short range.
165 * Total virial is computed in global_stat, called from do_md
167 const rvec* f = as_rvec_array(forceWithShiftForces.force().data());
168 f_calc_vir(start, start + homenr, x, f, vir_part, box);
169 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
173 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
177 static void pull_potential_wrapper(const t_commrec* cr,
178 const t_inputrec* ir,
180 gmx::ArrayRef<const gmx::RVec> x,
181 gmx::ForceWithVirial* force,
182 const t_mdatoms* mdatoms,
183 gmx_enerdata_t* enerd,
187 gmx_wallcycle_t wcycle)
192 /* Calculate the center of mass forces, this requires communication,
193 * which is why pull_potential is called close to other communication.
195 wallcycle_start(wcycle, ewcPULLPOT);
196 set_pbc(&pbc, ir->pbcType, box);
198 enerd->term[F_COM_PULL] += pull_potential(pull_work, mdatoms, &pbc, cr, t, lambda[efptRESTRAINT],
199 as_rvec_array(x.data()), force, &dvdl);
200 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
201 for (auto& dhdl : enerd->dhdlLambda)
205 wallcycle_stop(wcycle, ewcPULLPOT);
208 static void pme_receive_force_ener(t_forcerec* fr,
210 gmx::ForceWithVirial* forceWithVirial,
211 gmx_enerdata_t* enerd,
212 bool useGpuPmePpComms,
213 bool receivePmeForceToGpu,
214 gmx_wallcycle_t wcycle)
216 real e_q, e_lj, dvdl_q, dvdl_lj;
217 float cycles_ppdpme, cycles_seppme;
219 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
220 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
222 /* In case of node-splitting, the PP nodes receive the long-range
223 * forces, virial and energy from the PME nodes here.
225 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
228 gmx_pme_receive_f(fr->pmePpCommGpu.get(), cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
229 useGpuPmePpComms, receivePmeForceToGpu, &cycles_seppme);
230 enerd->term[F_COUL_RECIP] += e_q;
231 enerd->term[F_LJ_RECIP] += e_lj;
232 enerd->dvdl_lin[efptCOUL] += dvdl_q;
233 enerd->dvdl_lin[efptVDW] += dvdl_lj;
235 for (auto& dhdl : enerd->dhdlLambda)
237 dhdl += dvdl_q + dvdl_lj;
242 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
244 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
247 static void print_large_forces(FILE* fp,
252 ArrayRef<const RVec> x,
255 real force2Tolerance = gmx::square(forceTolerance);
256 gmx::index numNonFinite = 0;
257 for (int i = 0; i < md->homenr; i++)
259 real force2 = norm2(f[i]);
260 bool nonFinite = !std::isfinite(force2);
261 if (force2 >= force2Tolerance || nonFinite)
263 fprintf(fp, "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n", step,
264 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
271 if (numNonFinite > 0)
273 /* Note that with MPI this fatal call on one rank might interrupt
274 * the printing on other ranks. But we can only avoid that with
275 * an expensive MPI barrier that we would need at each step.
277 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
281 static void post_process_forces(const t_commrec* cr,
284 gmx_wallcycle_t wcycle,
286 ArrayRef<const RVec> x,
287 ForceOutputs* forceOutputs,
289 const t_mdatoms* mdatoms,
290 const t_forcerec* fr,
291 gmx::VirtualSitesHandler* vsite,
292 const StepWorkload& stepWork)
294 rvec* f = as_rvec_array(forceOutputs->forceWithShiftForces().force().data());
296 if (forceOutputs->haveForceWithVirial())
298 auto& forceWithVirial = forceOutputs->forceWithVirial();
302 /* Spread the mesh force on virtual sites to the other particles...
303 * This is parallellized. MPI communication is performed
304 * if the constructing atoms aren't local.
306 const gmx::VirtualSitesHandler::VirialHandling virialHandling =
307 (stepWork.computeVirial ? gmx::VirtualSitesHandler::VirialHandling::NonLinear
308 : gmx::VirtualSitesHandler::VirialHandling::None);
309 matrix virial = { { 0 } };
310 vsite->spreadForces(x, forceWithVirial.force_, virialHandling, {}, virial, nrnb, box, wcycle);
311 forceWithVirial.addVirialContribution(virial);
314 if (stepWork.computeVirial)
316 /* Now add the forces, this is local */
317 sum_forces(f, forceWithVirial.force_);
319 /* Add the direct virial contributions */
321 forceWithVirial.computeVirial_,
322 "forceWithVirial should request virial computation when we request the virial");
323 m_add(vir_force, forceWithVirial.getVirial(), vir_force);
327 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
332 if (fr->print_force >= 0)
334 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
338 static void do_nb_verlet(t_forcerec* fr,
339 const interaction_const_t* ic,
340 gmx_enerdata_t* enerd,
341 const StepWorkload& stepWork,
342 const InteractionLocality ilocality,
346 gmx_wallcycle_t wcycle)
348 if (!stepWork.computeNonbondedForces)
350 /* skip non-bonded calculation */
354 nonbonded_verlet_t* nbv = fr->nbv.get();
356 /* GPU kernel launch overhead is already timed separately */
357 if (fr->cutoff_scheme != ecutsVERLET)
359 gmx_incons("Invalid cut-off scheme passed!");
364 /* When dynamic pair-list pruning is requested, we need to prune
365 * at nstlistPrune steps.
367 if (nbv->isDynamicPruningStepCpu(step))
369 /* Prune the pair-list beyond fr->ic->rlistPrune using
370 * the current coordinates of the atoms.
372 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
373 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
374 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
378 nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
381 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
383 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, v.ssize());
385 /* Note that we would like to avoid this conditional by putting it
386 * into the omp pragma instead, but then we still take the full
387 * omp parallel for overhead (at least with gcc5).
389 if (!useOpenmpThreading || nth == 1)
398 #pragma omp parallel for num_threads(nth) schedule(static)
399 for (gmx::index i = 0; i < v.ssize(); i++)
406 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
408 * \param groupOptions Group options, containing T-coupling options
410 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
412 real nrdfCoupled = 0;
413 real nrdfUncoupled = 0;
414 real kineticEnergy = 0;
415 for (int g = 0; g < groupOptions.ngtc; g++)
417 if (groupOptions.tau_t[g] >= 0)
419 nrdfCoupled += groupOptions.nrdf[g];
420 kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * BOLTZ;
424 nrdfUncoupled += groupOptions.nrdf[g];
428 /* This conditional with > also catches nrdf=0 */
429 if (nrdfCoupled > nrdfUncoupled)
431 return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
439 /*! \brief This routine checks that the potential energy is finite.
441 * Always checks that the potential energy is finite. If step equals
442 * inputrec.init_step also checks that the magnitude of the potential energy
443 * is reasonable. Terminates with a fatal error when a check fails.
444 * Note that passing this check does not guarantee finite forces,
445 * since those use slightly different arithmetics. But in most cases
446 * there is just a narrow coordinate range where forces are not finite
447 * and energies are finite.
449 * \param[in] step The step number, used for checking and printing
450 * \param[in] enerd The energy data; the non-bonded group energies need to be added to
451 * enerd.term[F_EPOT] before calling this routine \param[in] inputrec The input record
453 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
455 /* Threshold valid for comparing absolute potential energy against
456 * the kinetic energy. Normally one should not consider absolute
457 * potential energy values, but with a factor of one million
458 * we should never get false positives.
460 constexpr real c_thresholdFactor = 1e6;
462 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
463 real averageKineticEnergy = 0;
464 /* We only check for large potential energy at the initial step,
465 * because that is by far the most likely step for this too occur
466 * and because computing the average kinetic energy is not free.
467 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
468 * before they become NaN.
470 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
472 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
475 if (energyIsNotFinite
476 || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
481 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
482 "contributions to the energy are %g and %g, respectively. A %s potential energy "
483 "can be caused by overlapping interactions in bonded interactions or very large%s "
484 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
485 "configuration, incorrect interactions or parameters in the topology.",
486 step, enerd.term[F_EPOT], energyIsNotFinite ? "not finite" : "extremely high",
487 enerd.term[F_LJ], enerd.term[F_COUL_SR],
488 energyIsNotFinite ? "non-finite" : "very high", energyIsNotFinite ? " or Nan" : "");
492 /*! \brief Return true if there are special forces computed this step.
494 * The conditionals exactly correspond to those in computeSpecialForces().
496 static bool haveSpecialForces(const t_inputrec& inputrec,
497 const gmx::ForceProviders& forceProviders,
498 const pull_t* pull_work,
499 const bool computeForces,
503 return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
504 (inputrec.bPull && pull_have_potential(pull_work)) || // pull
505 inputrec.bRot || // enforced rotation
506 (ed != nullptr) || // flooding
507 (inputrec.bIMD && computeForces)); // IMD
510 /*! \brief Compute forces and/or energies for special algorithms
512 * The intention is to collect all calls to algorithms that compute
513 * forces on local atoms only and that do not contribute to the local
514 * virial sum (but add their virial contribution separately).
515 * Eventually these should likely all become ForceProviders.
516 * Within this function the intention is to have algorithms that do
517 * global communication at the end, so global barriers within the MD loop
518 * are as close together as possible.
520 * \param[in] fplog The log file
521 * \param[in] cr The communication record
522 * \param[in] inputrec The input record
523 * \param[in] awh The Awh module (nullptr if none in use).
524 * \param[in] enforcedRotation Enforced rotation module.
525 * \param[in] imdSession The IMD session
526 * \param[in] pull_work The pull work structure.
527 * \param[in] step The current MD step
528 * \param[in] t The current time
529 * \param[in,out] wcycle Wallcycle accounting struct
530 * \param[in,out] forceProviders Pointer to a list of force providers
531 * \param[in] box The unit cell
532 * \param[in] x The coordinates
533 * \param[in] mdatoms Per atom properties
534 * \param[in] lambda Array of free-energy lambda values
535 * \param[in] stepWork Step schedule flags
536 * \param[in,out] forceWithVirial Force and virial buffers
537 * \param[in,out] enerd Energy buffer
538 * \param[in,out] ed Essential dynamics pointer
539 * \param[in] didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
541 * \todo Remove didNeighborSearch, which is used incorrectly.
542 * \todo Convert all other algorithms called here to ForceProviders.
544 static void computeSpecialForces(FILE* fplog,
546 const t_inputrec* inputrec,
548 gmx_enfrot* enforcedRotation,
549 gmx::ImdSession* imdSession,
553 gmx_wallcycle_t wcycle,
554 gmx::ForceProviders* forceProviders,
556 gmx::ArrayRef<const gmx::RVec> x,
557 const t_mdatoms* mdatoms,
559 const StepWorkload& stepWork,
560 gmx::ForceWithVirial* forceWithVirial,
561 gmx_enerdata_t* enerd,
563 bool didNeighborSearch)
565 /* NOTE: Currently all ForceProviders only provide forces.
566 * When they also provide energies, remove this conditional.
568 if (stepWork.computeForces)
570 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
571 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
573 /* Collect forces from modules */
574 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
577 if (inputrec->bPull && pull_have_potential(pull_work))
579 pull_potential_wrapper(cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work,
584 enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(
585 inputrec->pbcType, *mdatoms, box, forceWithVirial, t, step, wcycle, fplog);
589 rvec* f = as_rvec_array(forceWithVirial->force_.data());
591 /* Add the forces from enforced rotation potentials (if any) */
594 wallcycle_start(wcycle, ewcROTadd);
595 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
596 wallcycle_stop(wcycle, ewcROTadd);
601 /* Note that since init_edsam() is called after the initialization
602 * of forcerec, edsam doesn't request the noVirSum force buffer.
603 * Thus if no other algorithm (e.g. PME) requires it, the forces
604 * here will contribute to the virial.
606 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
609 /* Add forces from interactive molecular dynamics (IMD), if any */
610 if (inputrec->bIMD && stepWork.computeForces)
612 imdSession->applyForces(f);
616 /*! \brief Launch the prepare_step and spread stages of PME GPU.
618 * \param[in] pmedata The PME structure
619 * \param[in] box The box matrix
620 * \param[in] stepWork Step schedule flags
621 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in
622 * the device memory. \param[in] wcycle The wallcycle structure
624 static inline void launchPmeGpuSpread(gmx_pme_t* pmedata,
626 const StepWorkload& stepWork,
627 GpuEventSynchronizer* xReadyOnDevice,
628 gmx_wallcycle_t wcycle)
630 pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
631 pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle);
634 /*! \brief Launch the FFT and gather stages of PME GPU
636 * This function only implements setting the output forces (no accumulation).
638 * \param[in] pmedata The PME structure
639 * \param[in] wcycle The wallcycle structure
640 * \param[in] stepWork Step schedule flags
642 static void launchPmeGpuFftAndGather(gmx_pme_t* pmedata, gmx_wallcycle_t wcycle, const gmx::StepWorkload& stepWork)
644 pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
645 pme_gpu_launch_gather(pmedata, wcycle);
649 * Polling wait for either of the PME or nonbonded GPU tasks.
651 * Instead of a static order in waiting for GPU tasks, this function
652 * polls checking which of the two tasks completes first, and does the
653 * associated force buffer reduction overlapped with the other task.
654 * By doing that, unlike static scheduling order, it can always overlap
655 * one of the reductions, regardless of the GPU task completion order.
657 * \param[in] nbv Nonbonded verlet structure
658 * \param[in,out] pmedata PME module data
659 * \param[in,out] forceOutputs Output buffer for the forces and virial
660 * \param[in,out] enerd Energy data structure results are reduced into
661 * \param[in] stepWork Step schedule flags
662 * \param[in] wcycle The wallcycle structure
664 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
666 gmx::ForceOutputs* forceOutputs,
667 gmx_enerdata_t* enerd,
668 const StepWorkload& stepWork,
669 gmx_wallcycle_t wcycle)
671 bool isPmeGpuDone = false;
672 bool isNbGpuDone = false;
675 gmx::ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
676 gmx::ForceWithVirial& forceWithVirial = forceOutputs->forceWithVirial();
678 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
680 while (!isPmeGpuDone || !isNbGpuDone)
684 GpuTaskCompletion completionType =
685 (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
686 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, stepWork, wcycle, &forceWithVirial,
687 enerd, completionType);
692 GpuTaskCompletion completionType =
693 (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
694 isNbGpuDone = Nbnxm::gpu_try_finish_task(
695 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
696 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(),
697 completionType, wcycle);
701 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShiftForces.force());
707 /*! \brief Set up the different force buffers; also does clearing.
709 * \param[in] forceHelperBuffers Helper force buffers
710 * \param[in] pull_work The pull work object.
711 * \param[in] inputrec input record
712 * \param[in] force force array
713 * \param[in] stepWork Step schedule flags
714 * \param[out] wcycle wallcycle recording structure
716 * \returns Cleared force output structure
718 static ForceOutputs setupForceOutputs(ForceHelperBuffers* forceHelperBuffers,
720 const t_inputrec& inputrec,
721 gmx::ArrayRefWithPadding<gmx::RVec> force,
722 const StepWorkload& stepWork,
723 gmx_wallcycle_t wcycle)
725 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
727 /* NOTE: We assume fr->shiftForces is all zeros here */
728 gmx::ForceWithShiftForces forceWithShiftForces(force, stepWork.computeVirial,
729 forceHelperBuffers->shiftForces());
731 if (stepWork.computeForces)
733 /* Clear the short- and long-range forces */
734 clearRVecs(forceWithShiftForces.force(), true);
736 /* Clear the shift forces */
737 clearRVecs(forceWithShiftForces.shiftForces(), false);
740 /* If we need to compute the virial, we might need a separate
741 * force buffer for algorithms for which the virial is calculated
742 * directly, such as PME. Otherwise, forceWithVirial uses the
743 * the same force (f in legacy calls) buffer as other algorithms.
745 const bool useSeparateForceWithVirialBuffer =
746 (stepWork.computeForces
747 && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
748 /* forceWithVirial uses the local atom range only */
749 gmx::ForceWithVirial forceWithVirial(
750 useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
751 : force.unpaddedArrayRef(),
752 stepWork.computeVirial);
754 if (useSeparateForceWithVirialBuffer)
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.
761 clearRVecs(forceWithVirial.force_, true);
764 if (inputrec.bPull && pull_have_constraint(pull_work))
766 clear_pull_forces(pull_work);
769 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
771 return ForceOutputs(forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(),
776 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
778 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec& inputrec,
779 const t_forcerec& fr,
780 const pull_t* pull_work,
782 const InteractionDefinitions& idef,
784 const t_mdatoms& mdatoms,
785 const SimulationWorkload& simulationWork,
786 const StepWorkload& stepWork)
788 DomainLifetimeWorkload domainWork;
789 // Note that haveSpecialForces is constant over the whole run
790 domainWork.haveSpecialForces =
791 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
792 domainWork.haveCpuBondedWork = haveCpuBondeds(fr);
793 domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
794 domainWork.haveRestraintsWork = haveRestraints(idef, fcd);
795 domainWork.haveCpuListedForceWork = haveCpuListedForces(fr, idef, fcd);
796 // Note that haveFreeEnergyWork is constant over the whole run
797 domainWork.haveFreeEnergyWork = (fr.efep != efepNO && mdatoms.nPerturbed != 0);
798 // We assume we have local force work if there are CPU
799 // force tasks including PME or nonbondeds.
800 domainWork.haveCpuLocalForceWork =
801 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
802 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
803 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
808 /*! \brief Set up force flag stuct from the force bitmask.
810 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
811 * \param[in] isNonbondedOn Global override, if false forces to turn off all nonbonded calculation.
812 * \param[in] simulationWork Simulation workload description.
813 * \param[in] rankHasPmeDuty If this rank computes PME.
815 * \returns New Stepworkload description.
817 static StepWorkload setupStepWorkload(const int legacyFlags,
818 const bool isNonbondedOn,
819 const SimulationWorkload& simulationWork,
820 const bool rankHasPmeDuty)
823 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
824 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
825 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
826 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
827 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
828 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
829 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
830 flags.computeNonbondedForces = ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && isNonbondedOn;
831 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
833 if (simulationWork.useGpuBufferOps)
835 GMX_ASSERT(simulationWork.useGpuNonbonded,
836 "Can only offload buffer ops if nonbonded computation is also offloaded");
838 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
839 // on virial steps the CPU reduction path is taken
840 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
841 flags.useGpuPmeFReduction = flags.useGpuFBufferOps
842 && (simulationWork.useGpuPme
843 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication));
849 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
851 * TODO: eliminate \p useGpuPmeOnThisRank when this is
852 * incorporated in DomainLifetimeWorkload.
854 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
855 gmx::GpuBonded* gpuBonded,
857 gmx_enerdata_t* enerd,
858 const gmx::MdrunScheduleWorkload& runScheduleWork,
859 bool useGpuPmeOnThisRank,
861 gmx_wallcycle_t wcycle)
863 if (runScheduleWork.simulationWork.useGpuNonbonded)
865 /* Launch pruning before buffer clearing because the API overhead of the
866 * clear kernel launches can leave the GPU idle while it could be running
869 if (nbv->isDynamicPruningStepGpu(step))
871 nbv->dispatchPruneKernelGpu(step);
874 /* now clear the GPU outputs while we finish the step on the CPU */
875 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
876 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
877 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
878 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
879 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
882 if (useGpuPmeOnThisRank)
884 pme_gpu_reinit_computation(pmedata, wcycle);
887 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
889 // in principle this should be included in the DD balancing region,
890 // but generally it is infrequent so we'll omit it for the sake of
892 gpuBonded->waitAccumulateEnergyTerms(enerd);
894 gpuBonded->clearEnergies();
898 //! \brief Data structure to hold dipole-related data and staging arrays
901 //! Dipole staging for fast summing over MPI
902 gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
903 //! Dipole staging for states A and B (index 0 and 1 resp.)
904 gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
908 static void reduceAndUpdateMuTot(DipoleData* dipoleData,
910 const bool haveFreeEnergy,
911 gmx::ArrayRef<const real> lambda,
913 const DDBalanceRegionHandler& ddBalanceRegionHandler)
917 gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
918 ddBalanceRegionHandler.reopenRegionCpu();
920 for (int i = 0; i < 2; i++)
922 for (int j = 0; j < DIM; j++)
924 dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
930 copy_rvec(dipoleData->muStateAB[0], muTotal);
934 for (int j = 0; j < DIM; j++)
936 muTotal[j] = (1.0 - lambda[efptCOUL]) * dipoleData->muStateAB[0][j]
937 + lambda[efptCOUL] * dipoleData->muStateAB[1][j];
942 void do_force(FILE* fplog,
944 const gmx_multisim_t* ms,
945 const t_inputrec* inputrec,
947 gmx_enfrot* enforcedRotation,
948 gmx::ImdSession* imdSession,
952 gmx_wallcycle_t wcycle,
953 const gmx_localtop_t* top,
955 gmx::ArrayRefWithPadding<gmx::RVec> x,
957 gmx::ArrayRefWithPadding<gmx::RVec> force,
959 const t_mdatoms* mdatoms,
960 gmx_enerdata_t* enerd,
962 gmx::ArrayRef<real> lambda,
964 gmx::MdrunScheduleWorkload* runScheduleWork,
965 gmx::VirtualSitesHandler* vsite,
970 const DDBalanceRegionHandler& ddBalanceRegionHandler)
972 GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
973 "The size of the force buffer should be at least the number of atoms to compute "
976 nonbonded_verlet_t* nbv = fr->nbv.get();
977 interaction_const_t* ic = fr->ic;
978 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
980 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
983 runScheduleWork->stepWork = setupStepWorkload(legacyFlags, fr->bNonbonded, simulationWork,
984 thisRankHasDuty(cr, DUTY_PME));
985 const StepWorkload& stepWork = runScheduleWork->stepWork;
988 const bool useGpuPmeOnThisRank = simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME);
990 /* At a search step we need to start the first balancing region
991 * somewhere early inside the step after communication during domain
992 * decomposition (and not during the previous step as usual).
994 if (stepWork.doNeighborSearch)
996 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
999 clear_mat(vir_force);
1001 if (fr->pbcType != PbcType::No)
1003 /* Compute shift vectors every step,
1004 * because of pressure coupling or box deformation!
1006 if (stepWork.haveDynamicBox && stepWork.stateChanged)
1008 calc_shifts(box, fr->shift_vec);
1011 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1012 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1015 put_atoms_in_box_omp(fr->pbcType, box, x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1016 gmx_omp_nthreads_get(emntDefault));
1017 inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1021 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1023 const bool pmeSendCoordinatesFromGpu =
1024 GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1025 const bool reinitGpuPmePpComms =
1026 GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1028 const auto localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1029 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1030 AtomLocality::Local, simulationWork, stepWork)
1033 // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1034 // Otherwise the send will occur after H2D coordinate transfer.
1035 if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu)
1037 /* Send particle coordinates to the pme nodes */
1038 if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1040 GMX_RELEASE_ASSERT(false,
1041 "GPU update and separate PME ranks are only supported with GPU "
1042 "direct communication!");
1043 // TODO: when this code-path becomes supported add:
1044 // stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1047 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1048 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1049 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1050 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1053 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1054 // The local coordinates can be copied right away.
1055 // NOTE: Consider moving this copy to right after they are updated and constrained,
1056 // if the later is not offloaded.
1057 if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1059 if (stepWork.doNeighborSearch)
1061 // TODO refactor this to do_md, after partitioning.
1062 stateGpu->reinit(mdatoms->homenr,
1063 cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1064 if (useGpuPmeOnThisRank)
1066 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1067 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1070 // We need to copy coordinates when:
1071 // 1. Update is not offloaded
1072 // 2. The buffers were reinitialized on search step
1073 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1075 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1076 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1080 // TODO Update this comment when introducing SimulationWorkload
1082 // The conditions for gpuHaloExchange e.g. using GPU buffer
1083 // operations were checked before construction, so here we can
1084 // just use it and assert upon any conditions.
1085 const bool ddUsesGpuDirectCommunication =
1086 ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange.empty()));
1087 GMX_ASSERT(!ddUsesGpuDirectCommunication || stepWork.useGpuXBufferOps,
1088 "Must use coordinate buffer ops with GPU halo exchange");
1089 const bool useGpuForcesHaloExchange = ddUsesGpuDirectCommunication && stepWork.useGpuFBufferOps;
1091 // Copy coordinate from the GPU if update is on the GPU and there
1092 // are forces to be computed on the CPU, or for the computation of
1093 // virial, or if host-side data will be transferred from this task
1094 // to a remote task for halo exchange or PME-PP communication. At
1095 // search steps the current coordinates are already on the host,
1096 // hence copy is not needed.
1097 const bool haveHostPmePpComms =
1098 !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1099 const bool haveHostHaloExchangeComms = havePPDomainDecomposition(cr) && !ddUsesGpuDirectCommunication;
1101 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1102 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1103 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1104 || haveHostPmePpComms || haveHostHaloExchangeComms))
1106 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1107 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1108 haveCopiedXFromGpu = true;
1111 // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1112 // Otherwise the send will occur before the H2D coordinate transfer.
1113 if (!thisRankHasDuty(cr, DUTY_PME) && pmeSendCoordinatesFromGpu)
1115 /* Send particle coordinates to the pme nodes */
1116 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1117 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1118 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1119 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1122 if (useGpuPmeOnThisRank)
1124 launchPmeGpuSpread(fr->pmedata, box, stepWork, localXReadyOnDevice, wcycle);
1127 /* do gridding for pair search */
1128 if (stepWork.doNeighborSearch)
1130 if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1132 fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1136 // - vzero is constant, do we need to pass it?
1137 // - box_diag should be passed directly to nbnxn_put_on_grid
1143 box_diag[XX] = box[XX][XX];
1144 box_diag[YY] = box[YY][YY];
1145 box_diag[ZZ] = box[ZZ][ZZ];
1147 wallcycle_start(wcycle, ewcNS);
1148 if (!DOMAINDECOMP(cr))
1150 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1151 nbnxn_put_on_grid(nbv, box, 0, vzero, box_diag, nullptr, { 0, mdatoms->homenr }, -1,
1152 fr->cginfo, x.unpaddedArrayRef(), 0, nullptr);
1153 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1157 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1158 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1159 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1162 nbv->setAtomProperties(*mdatoms, fr->cginfo);
1164 wallcycle_stop(wcycle, ewcNS);
1166 /* initialize the GPU nbnxm atom data and bonded data structures */
1167 if (simulationWork.useGpuNonbonded)
1169 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1171 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1172 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1173 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1177 /* Now we put all atoms on the grid, we can assign bonded
1178 * interactions to the GPU, where the grid order is
1179 * needed. Also the xq, f and fshift device buffers have
1180 * been reallocated if needed, so the bonded code can
1181 * learn about them. */
1182 // TODO the xq, f, and fshift buffers are now shared
1183 // resources, so they should be maintained by a
1184 // higher-level object than the nb module.
1185 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(
1186 nbv->getGridIndices(), top->idef, Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1187 Nbnxm::gpu_get_f(nbv->gpu_nbv), Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1189 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1192 // Need to run after the GPU-offload bonded interaction lists
1193 // are set up to be able to determine whether there is bonded work.
1194 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1195 *inputrec, *fr, pull_work, ed, top->idef, *fcd, *mdatoms, simulationWork, stepWork);
1197 wallcycle_start_nocount(wcycle, ewcNS);
1198 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1199 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1200 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1202 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1204 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1205 wallcycle_stop(wcycle, ewcNS);
1207 if (stepWork.useGpuXBufferOps)
1209 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1211 // For force buffer ops, we use the below conditon rather than
1212 // useGpuFBufferOps to ensure that init is performed even if this
1213 // NS step is also a virial step (on which f buf ops are deactivated).
1214 if (simulationWork.useGpuBufferOps && simulationWork.useGpuNonbonded && (GMX_GPU == GMX_GPU_CUDA))
1216 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1217 nbv->atomdata_init_add_nbat_f_to_f_gpu(stateGpu->fReducedOnDevice());
1220 else if (!EI_TPI(inputrec->eI))
1222 if (stepWork.useGpuXBufferOps)
1224 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1225 nbv->convertCoordinatesGpu(AtomLocality::Local, false, stateGpu->getCoordinates(),
1226 localXReadyOnDevice);
1230 if (simulationWork.useGpuUpdate)
1232 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1233 GMX_ASSERT(haveCopiedXFromGpu,
1234 "a wait should only be triggered if copy has been scheduled");
1235 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1237 nbv->convertCoordinates(AtomLocality::Local, false, x.unpaddedArrayRef());
1241 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1243 if (simulationWork.useGpuNonbonded)
1245 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1247 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1249 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1250 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1251 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1253 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1255 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1256 // with X buffer ops offloaded to the GPU on all but the search steps
1258 // bonded work not split into separate local and non-local, so with DD
1259 // we can only launch the kernel after non-local coordinates have been received.
1260 if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1262 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1263 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1264 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1267 /* launch local nonbonded work on GPU */
1268 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1269 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1270 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1271 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1274 if (useGpuPmeOnThisRank)
1276 // In PME GPU and mixed mode we launch FFT / gather after the
1277 // X copy/transform to allow overlap as well as after the GPU NB
1278 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1279 // the nonbonded kernel.
1280 launchPmeGpuFftAndGather(fr->pmedata, wcycle, stepWork);
1283 /* Communicate coordinates and sum dipole if necessary +
1284 do non-local pair search */
1285 if (havePPDomainDecomposition(cr))
1287 if (stepWork.doNeighborSearch)
1289 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1290 wallcycle_start_nocount(wcycle, ewcNS);
1291 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1292 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1293 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1295 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1296 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1297 wallcycle_stop(wcycle, ewcNS);
1298 // TODO refactor this GPU halo exchange re-initialisation
1299 // to location in do_md where GPU halo exchange is
1300 // constructed at partitioning, after above stateGpu
1301 // re-initialization has similarly been refactored
1302 if (ddUsesGpuDirectCommunication)
1304 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1309 if (ddUsesGpuDirectCommunication)
1311 // The following must be called after local setCoordinates (which records an event
1312 // when the coordinate data has been copied to the device).
1313 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1315 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1317 // non-local part of coordinate buffer must be copied back to host for CPU work
1318 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1323 // Note: GPU update + DD without direct communication is not supported,
1324 // a waitCoordinatesReadyOnHost() should be issued if it will be.
1325 GMX_ASSERT(!simulationWork.useGpuUpdate,
1326 "GPU update is not supported with CPU halo exchange");
1327 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1330 if (stepWork.useGpuXBufferOps)
1332 if (!useGpuPmeOnThisRank && !ddUsesGpuDirectCommunication)
1334 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1336 nbv->convertCoordinatesGpu(AtomLocality::NonLocal, false, stateGpu->getCoordinates(),
1337 stateGpu->getCoordinatesReadyOnDeviceEvent(
1338 AtomLocality::NonLocal, simulationWork, stepWork));
1342 nbv->convertCoordinates(AtomLocality::NonLocal, false, x.unpaddedArrayRef());
1346 if (simulationWork.useGpuNonbonded)
1348 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1350 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1352 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1353 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1354 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1357 if (domainWork.haveGpuBondedWork)
1359 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1360 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1361 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1364 /* launch non-local nonbonded tasks on GPU */
1365 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1366 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1368 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1370 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1374 if (simulationWork.useGpuNonbonded)
1376 /* launch D2H copy-back F */
1377 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1378 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1380 if (havePPDomainDecomposition(cr))
1382 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1384 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1385 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1387 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1389 fr->gpuBonded->launchEnergyTransfer();
1391 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1394 gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1395 if (fr->wholeMoleculeTransform)
1397 xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1400 DipoleData dipoleData;
1402 if (simulationWork.computeMuTot)
1404 const int start = 0;
1406 /* Calculate total (local) dipole moment in a temporary common array.
1407 * This makes it possible to sum them over nodes faster.
1409 gmx::ArrayRef<const gmx::RVec> xRef =
1410 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1411 calc_mu(start, mdatoms->homenr, xRef, mdatoms->chargeA, mdatoms->chargeB,
1412 mdatoms->nChargePerturbed, dipoleData.muStaging[0], dipoleData.muStaging[1]);
1414 reduceAndUpdateMuTot(&dipoleData, cr, (fr->efep != efepNO), lambda, muTotal, ddBalanceRegionHandler);
1417 /* Reset energies */
1418 reset_enerdata(enerd);
1420 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1422 wallcycle_start(wcycle, ewcPPDURINGPME);
1423 dd_force_flop_start(cr->dd, nrnb);
1426 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1427 // this wait ensures that the D2H transfer is complete.
1428 if ((simulationWork.useGpuUpdate)
1429 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1431 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1436 wallcycle_start(wcycle, ewcROT);
1437 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step,
1438 stepWork.doNeighborSearch);
1439 wallcycle_stop(wcycle, ewcROT);
1442 /* Start the force cycle counter.
1443 * Note that a different counter is used for dynamic load balancing.
1445 wallcycle_start(wcycle, ewcFORCE);
1447 // Set up and clear force outputs.
1448 // We use std::move to keep the compiler happy, it has no effect.
1449 ForceOutputs forceOut = setupForceOutputs(fr->forceHelperBuffers.get(), pull_work, *inputrec,
1450 std::move(force), stepWork, wcycle);
1452 /* We calculate the non-bonded forces, when done on the CPU, here.
1453 * We do this before calling do_force_lowlevel, because in that
1454 * function, the listed forces are calculated before PME, which
1455 * does communication. With this order, non-bonded and listed
1456 * force calculation imbalance can be balanced out by the domain
1457 * decomposition load balancing.
1460 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1462 if (!useOrEmulateGpuNb)
1464 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1467 if (fr->efep != efepNO)
1469 /* Calculate the local and non-local free energy interactions here.
1470 * Happens here on the CPU both with and without GPU.
1472 nbv->dispatchFreeEnergyKernel(InteractionLocality::Local, fr,
1473 as_rvec_array(x.unpaddedArrayRef().data()),
1474 &forceOut.forceWithShiftForces(), *mdatoms, inputrec->fepvals,
1475 lambda.data(), enerd, stepWork, nrnb);
1477 if (havePPDomainDecomposition(cr))
1479 nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal, fr,
1480 as_rvec_array(x.unpaddedArrayRef().data()),
1481 &forceOut.forceWithShiftForces(), *mdatoms,
1482 inputrec->fepvals, lambda.data(), enerd, stepWork, nrnb);
1486 if (!useOrEmulateGpuNb)
1488 if (havePPDomainDecomposition(cr))
1490 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1494 if (stepWork.computeForces)
1496 /* Add all the non-bonded force to the normal force array.
1497 * This can be split into a local and a non-local part when overlapping
1498 * communication with calculation with domain decomposition.
1500 wallcycle_stop(wcycle, ewcFORCE);
1501 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All, forceOut.forceWithShiftForces().force());
1502 wallcycle_start_nocount(wcycle, ewcFORCE);
1505 /* If there are multiple fshift output buffers we need to reduce them */
1506 if (stepWork.computeVirial)
1508 /* This is not in a subcounter because it takes a
1509 negligible and constant-sized amount of time */
1510 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat,
1511 forceOut.forceWithShiftForces().shiftForces());
1515 // TODO Force flags should include haveFreeEnergyWork for this domain
1516 if (ddUsesGpuDirectCommunication && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1518 /* Wait for non-local coordinate data to be copied from device */
1519 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1521 /* Compute the bonded and non-bonded energies and optionally forces */
1522 do_force_lowlevel(fr, inputrec, top->idef, cr, ms, nrnb, wcycle, mdatoms, x, xWholeMolecules,
1523 hist, &forceOut, enerd, fcd, box, lambda.data(),
1524 as_rvec_array(dipoleData.muStateAB), stepWork, ddBalanceRegionHandler);
1526 wallcycle_stop(wcycle, ewcFORCE);
1528 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation, imdSession, pull_work, step, t,
1529 wcycle, fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda.data(),
1530 stepWork, &forceOut.forceWithVirial(), enerd, ed, stepWork.doNeighborSearch);
1533 // Will store the amount of cycles spent waiting for the GPU that
1534 // will be later used in the DLB accounting.
1535 float cycles_wait_gpu = 0;
1536 if (useOrEmulateGpuNb)
1538 auto& forceWithShiftForces = forceOut.forceWithShiftForces();
1540 /* wait for non-local forces (or calculate in emulation mode) */
1541 if (havePPDomainDecomposition(cr))
1543 if (simulationWork.useGpuNonbonded)
1545 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1546 nbv->gpu_nbv, stepWork, AtomLocality::NonLocal, enerd->grpp.ener[egLJSR].data(),
1547 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(), wcycle);
1551 wallcycle_start_nocount(wcycle, ewcFORCE);
1552 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes,
1553 step, nrnb, wcycle);
1554 wallcycle_stop(wcycle, ewcFORCE);
1557 if (stepWork.useGpuFBufferOps)
1559 gmx::FixedCapacityVector<GpuEventSynchronizer*, 1> dependencyList;
1561 // TODO: move this into DomainLifetimeWorkload, including the second part of the
1562 // condition The bonded and free energy CPU tasks can have non-local force
1563 // contributions which are a dependency for the GPU force reduction.
1564 bool haveNonLocalForceContribInCpuBuffer =
1565 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1567 if (haveNonLocalForceContribInCpuBuffer)
1569 stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(),
1570 AtomLocality::NonLocal);
1571 dependencyList.push_back(stateGpu->getForcesReadyOnDeviceEvent(
1572 AtomLocality::NonLocal, stepWork.useGpuFBufferOps));
1575 nbv->atomdata_add_nbat_f_to_f_gpu(AtomLocality::NonLocal, stateGpu->getForces(),
1576 pme_gpu_get_device_f(fr->pmedata), dependencyList,
1577 false, haveNonLocalForceContribInCpuBuffer);
1578 if (!useGpuForcesHaloExchange)
1580 // copy from GPU input for dd_move_f()
1581 stateGpu->copyForcesFromGpu(forceOut.forceWithShiftForces().force(),
1582 AtomLocality::NonLocal);
1587 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
1591 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
1593 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
1598 if (havePPDomainDecomposition(cr))
1600 /* We are done with the CPU compute.
1601 * We will now communicate the non-local forces.
1602 * If we use a GPU this will overlap with GPU work, so in that case
1603 * we do not close the DD force balancing region here.
1605 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1607 if (stepWork.computeForces)
1610 if (useGpuForcesHaloExchange)
1612 if (domainWork.haveCpuLocalForceWork)
1614 stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(), AtomLocality::Local);
1616 communicateGpuHaloForces(*cr, domainWork.haveCpuLocalForceWork);
1620 if (stepWork.useGpuFBufferOps)
1622 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
1624 dd_move_f(cr->dd, &forceOut.forceWithShiftForces(), wcycle);
1629 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1630 // an alternating wait/reduction scheme.
1631 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
1632 && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
1633 if (alternateGpuWait)
1635 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &forceOut, enerd, stepWork, wcycle);
1638 if (!alternateGpuWait && useGpuPmeOnThisRank)
1640 pme_gpu_wait_and_reduce(fr->pmedata, stepWork, wcycle, &forceOut.forceWithVirial(), enerd);
1643 /* Wait for local GPU NB outputs on the non-alternating wait path */
1644 if (!alternateGpuWait && simulationWork.useGpuNonbonded)
1646 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1647 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1648 * but even with a step of 0.1 ms the difference is less than 1%
1651 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
1652 const float waitCycles = Nbnxm::gpu_wait_finish_task(
1653 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
1654 enerd->grpp.ener[egCOULSR].data(), forceOut.forceWithShiftForces().shiftForces(), wcycle);
1656 if (ddBalanceRegionHandler.useBalancingRegion())
1658 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1659 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
1661 /* We measured few cycles, it could be that the kernel
1662 * and transfer finished earlier and there was no actual
1663 * wait time, only API call overhead.
1664 * Then the actual time could be anywhere between 0 and
1665 * cycles_wait_est. We will use half of cycles_wait_est.
1667 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1669 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1673 if (fr->nbv->emulateGpu())
1675 // NOTE: emulation kernel is not included in the balancing region,
1676 // but emulation mode does not target performance anyway
1677 wallcycle_start_nocount(wcycle, ewcFORCE);
1678 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local,
1679 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes, step, nrnb, wcycle);
1680 wallcycle_stop(wcycle, ewcFORCE);
1683 // If on GPU PME-PP comms or GPU update path, receive forces from PME before GPU buffer ops
1684 // TODO refactor this and unify with below default-path call to the same function
1685 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME)
1686 && (simulationWork.useGpuPmePpCommunication || simulationWork.useGpuUpdate))
1688 /* In case of node-splitting, the PP nodes receive the long-range
1689 * forces, virial and energy from the PME nodes here.
1691 pme_receive_force_ener(fr, cr, &forceOut.forceWithVirial(), enerd,
1692 simulationWork.useGpuPmePpCommunication,
1693 stepWork.useGpuPmeFReduction, wcycle);
1697 /* Do the nonbonded GPU (or emulation) force buffer reduction
1698 * on the non-alternating path. */
1699 if (useOrEmulateGpuNb && !alternateGpuWait)
1701 // TODO simplify the below conditionals. Pass buffer and sync pointers at init stage rather than here. Unify getter fns for sameGPU/otherGPU cases.
1703 stepWork.useGpuPmeFReduction
1704 ? (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1705 : // PME force buffer on same GPU
1706 fr->pmePpCommGpu->getGpuForceStagingPtr()) // buffer received from other GPU
1707 : nullptr; // PME reduction not active on GPU
1709 GpuEventSynchronizer* const pmeSynchronizer =
1710 stepWork.useGpuPmeFReduction
1711 ? (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1712 : // PME force buffer on same GPU
1713 static_cast<GpuEventSynchronizer*>(
1714 fr->pmePpCommGpu->getForcesReadySynchronizer())) // buffer received from other GPU
1715 : nullptr; // PME reduction not active on GPU
1717 gmx::FixedCapacityVector<GpuEventSynchronizer*, 3> dependencyList;
1719 if (stepWork.useGpuPmeFReduction)
1721 dependencyList.push_back(pmeSynchronizer);
1724 gmx::ArrayRef<gmx::RVec> forceWithShift = forceOut.forceWithShiftForces().force();
1726 if (stepWork.useGpuFBufferOps)
1728 // Flag to specify whether the CPU force buffer has contributions to
1729 // local atoms. This depends on whether there are CPU-based force tasks
1730 // or when DD is active the halo exchange has resulted in contributions
1731 // from the non-local part.
1732 const bool haveLocalForceContribInCpuBuffer =
1733 (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
1735 // TODO: move these steps as early as possible:
1736 // - CPU f H2D should be as soon as all CPU-side forces are done
1737 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
1738 // before the next CPU task that consumes the forces: vsite spread or update)
1739 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
1740 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
1741 // These should be unified.
1742 if (haveLocalForceContribInCpuBuffer && !useGpuForcesHaloExchange)
1744 // Note: AtomLocality::All is used for the non-DD case because, as in this
1745 // case copyForcesToGpu() uses a separate stream, it allows overlap of
1746 // CPU force H2D with GPU force tasks on all streams including those in the
1747 // local stream which would otherwise be implicit dependencies for the
1748 // transfer and would not overlap.
1749 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
1751 stateGpu->copyForcesToGpu(forceWithShift, locality);
1752 dependencyList.push_back(
1753 stateGpu->getForcesReadyOnDeviceEvent(locality, stepWork.useGpuFBufferOps));
1755 if (useGpuForcesHaloExchange)
1757 dependencyList.push_back(cr->dd->gpuHaloExchange[0]->getForcesReadyOnDeviceEvent());
1759 nbv->atomdata_add_nbat_f_to_f_gpu(AtomLocality::Local, stateGpu->getForces(), pmeForcePtr,
1760 dependencyList, stepWork.useGpuPmeFReduction,
1761 haveLocalForceContribInCpuBuffer);
1762 // Copy forces to host if they are needed for update or if virtual sites are enabled.
1763 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
1764 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1765 // copy call done in sim_utils(...) for the output.
1766 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
1767 // they should not be copied in do_md(...) for the output.
1768 if (!simulationWork.useGpuUpdate || vsite)
1770 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
1771 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1776 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
1780 launchGpuEndOfStepTasks(nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork,
1781 useGpuPmeOnThisRank, step, wcycle);
1783 if (DOMAINDECOMP(cr))
1785 dd_force_flop_stop(cr->dd, nrnb);
1788 if (stepWork.computeForces)
1790 /* If we have NoVirSum forces, but we do not calculate the virial,
1791 * we sum fr->f_novirsum=forceOut.f later.
1793 if (vsite && !(fr->forceHelperBuffers->haveDirectVirialContributions() && !stepWork.computeVirial))
1795 auto f = forceOut.forceWithShiftForces().force();
1796 auto fshift = forceOut.forceWithShiftForces().shiftForces();
1797 const gmx::VirtualSitesHandler::VirialHandling virialHandling =
1798 (stepWork.computeVirial ? gmx::VirtualSitesHandler::VirialHandling::Pbc
1799 : gmx::VirtualSitesHandler::VirialHandling::None);
1800 vsite->spreadForces(x.unpaddedArrayRef(), f, virialHandling, fshift, nullptr, nrnb, box, wcycle);
1803 if (stepWork.computeVirial)
1805 /* Calculation of the virial must be done after vsites! */
1806 calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()),
1807 forceOut.forceWithShiftForces(), vir_force, box, nrnb, fr, inputrec->pbcType);
1811 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
1812 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
1813 && !simulationWork.useGpuUpdate)
1815 /* In case of node-splitting, the PP nodes receive the long-range
1816 * forces, virial and energy from the PME nodes here.
1818 pme_receive_force_ener(fr, cr, &forceOut.forceWithVirial(), enerd,
1819 simulationWork.useGpuPmePpCommunication, false, wcycle);
1822 if (stepWork.computeForces)
1824 post_process_forces(cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOut, vir_force,
1825 mdatoms, fr, vsite, stepWork);
1828 if (stepWork.computeEnergy)
1830 /* Sum the potential energy terms from group contributions */
1831 sum_epot(&(enerd->grpp), enerd->term);
1833 if (!EI_TPI(inputrec->eI))
1835 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1839 /* In case we don't have constraints and are using GPUs, the next balancing
1840 * region starts here.
1841 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1842 * virial calculation and COM pulling, is not thus not included in
1843 * the balance timing, which is ok as most tasks do communication.
1845 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);