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,2021, 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.
49 #include "gromacs/applied_forces/awh/awh.h"
50 #include "gromacs/domdec/dlbtiming.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/domdec/gpuhaloexchange.h"
54 #include "gromacs/domdec/partition.h"
55 #include "gromacs/essentialdynamics/edsam.h"
56 #include "gromacs/ewald/pme.h"
57 #include "gromacs/ewald/pme_pp.h"
58 #include "gromacs/ewald/pme_pp_comm_gpu.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
61 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
62 #include "gromacs/gmxlib/nrnb.h"
63 #include "gromacs/gpu_utils/gpu_utils.h"
64 #include "gromacs/imd/imd.h"
65 #include "gromacs/listed_forces/disre.h"
66 #include "gromacs/listed_forces/listed_forces_gpu.h"
67 #include "gromacs/listed_forces/listed_forces.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/dispersioncorrection.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/force.h"
80 #include "gromacs/mdlib/force_flags.h"
81 #include "gromacs/mdlib/forcerec.h"
82 #include "gromacs/mdlib/gmx_omp_nthreads.h"
83 #include "gromacs/mdlib/update.h"
84 #include "gromacs/mdlib/vsite.h"
85 #include "gromacs/mdlib/wall.h"
86 #include "gromacs/mdlib/wholemoleculetransform.h"
87 #include "gromacs/mdtypes/commrec.h"
88 #include "gromacs/mdtypes/enerdata.h"
89 #include "gromacs/mdtypes/forcebuffers.h"
90 #include "gromacs/mdtypes/forceoutput.h"
91 #include "gromacs/mdtypes/forcerec.h"
92 #include "gromacs/mdtypes/iforceprovider.h"
93 #include "gromacs/mdtypes/inputrec.h"
94 #include "gromacs/mdtypes/md_enums.h"
95 #include "gromacs/mdtypes/mdatom.h"
96 #include "gromacs/mdtypes/multipletimestepping.h"
97 #include "gromacs/mdtypes/simulation_workload.h"
98 #include "gromacs/mdtypes/state.h"
99 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
100 #include "gromacs/nbnxm/gpu_data_mgmt.h"
101 #include "gromacs/nbnxm/nbnxm.h"
102 #include "gromacs/nbnxm/nbnxm_gpu.h"
103 #include "gromacs/pbcutil/ishift.h"
104 #include "gromacs/pbcutil/pbc.h"
105 #include "gromacs/pulling/pull.h"
106 #include "gromacs/pulling/pull_rotation.h"
107 #include "gromacs/timing/cyclecounter.h"
108 #include "gromacs/timing/gpu_timing.h"
109 #include "gromacs/timing/wallcycle.h"
110 #include "gromacs/timing/wallcyclereporting.h"
111 #include "gromacs/timing/walltime_accounting.h"
112 #include "gromacs/topology/topology.h"
113 #include "gromacs/utility/arrayref.h"
114 #include "gromacs/utility/basedefinitions.h"
115 #include "gromacs/utility/cstringutil.h"
116 #include "gromacs/utility/exceptions.h"
117 #include "gromacs/utility/fatalerror.h"
118 #include "gromacs/utility/fixedcapacityvector.h"
119 #include "gromacs/utility/gmxassert.h"
120 #include "gromacs/utility/gmxmpi.h"
121 #include "gromacs/utility/logger.h"
122 #include "gromacs/utility/smalloc.h"
123 #include "gromacs/utility/strconvert.h"
124 #include "gromacs/utility/sysinfo.h"
126 #include "gpuforcereduction.h"
129 using gmx::AtomLocality;
130 using gmx::DomainLifetimeWorkload;
131 using gmx::ForceOutputs;
132 using gmx::ForceWithShiftForces;
133 using gmx::InteractionLocality;
135 using gmx::SimulationWorkload;
136 using gmx::StepWorkload;
138 // TODO: this environment variable allows us to verify before release
139 // that on less common architectures the total cost of polling is not larger than
140 // a blocking wait (so polling does not introduce overhead when the static
141 // PME-first ordering would suffice).
142 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
144 static void sum_forces(ArrayRef<RVec> f, ArrayRef<const RVec> forceToAdd)
146 GMX_ASSERT(f.size() >= forceToAdd.size(), "Accumulation buffer should be sufficiently large");
147 const int end = forceToAdd.size();
149 int gmx_unused nt = gmx_omp_nthreads_get(ModuleMultiThread::Default);
150 #pragma omp parallel for num_threads(nt) schedule(static)
151 for (int i = 0; i < end; i++)
153 rvec_inc(f[i], forceToAdd[i]);
157 static void calc_virial(int start,
160 const gmx::ForceWithShiftForces& forceWithShiftForces,
164 const t_forcerec* fr,
167 /* The short-range virial from surrounding boxes */
168 const rvec* fshift = as_rvec_array(forceWithShiftForces.shiftForces().data());
169 const rvec* shiftVecPointer = as_rvec_array(fr->shift_vec.data());
170 calc_vir(gmx::c_numShiftVectors, shiftVecPointer, fshift, vir_part, pbcType == PbcType::Screw, box);
171 inc_nrnb(nrnb, eNR_VIRIAL, gmx::c_numShiftVectors);
173 /* Calculate partial virial, for local atoms only, based on short range.
174 * Total virial is computed in global_stat, called from do_md
176 const rvec* f = as_rvec_array(forceWithShiftForces.force().data());
177 f_calc_vir(start, start + homenr, x, f, vir_part, box);
178 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
182 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
186 static void pull_potential_wrapper(const t_commrec* cr,
187 const t_inputrec& ir,
189 gmx::ArrayRef<const gmx::RVec> x,
190 gmx::ForceWithVirial* force,
191 const t_mdatoms* mdatoms,
192 gmx_enerdata_t* enerd,
196 gmx_wallcycle* wcycle)
201 /* Calculate the center of mass forces, this requires communication,
202 * which is why pull_potential is called close to other communication.
204 wallcycle_start(wcycle, WallCycleCounter::PullPot);
205 set_pbc(&pbc, ir.pbcType, box);
207 enerd->term[F_COM_PULL] +=
208 pull_potential(pull_work,
209 gmx::arrayRefFromArray(mdatoms->massT, mdatoms->nr),
213 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Restraint)],
217 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Restraint] += dvdl;
218 wallcycle_stop(wcycle, WallCycleCounter::PullPot);
221 static void pme_receive_force_ener(t_forcerec* fr,
223 gmx::ForceWithVirial* forceWithVirial,
224 gmx_enerdata_t* enerd,
225 bool useGpuPmePpComms,
226 bool receivePmeForceToGpu,
227 gmx_wallcycle* wcycle)
229 real e_q, e_lj, dvdl_q, dvdl_lj;
230 float cycles_ppdpme, cycles_seppme;
232 cycles_ppdpme = wallcycle_stop(wcycle, WallCycleCounter::PpDuringPme);
233 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
235 /* In case of node-splitting, the PP nodes receive the long-range
236 * forces, virial and energy from the PME nodes here.
238 wallcycle_start(wcycle, WallCycleCounter::PpPmeWaitRecvF);
241 gmx_pme_receive_f(fr->pmePpCommGpu.get(),
249 receivePmeForceToGpu,
251 enerd->term[F_COUL_RECIP] += e_q;
252 enerd->term[F_LJ_RECIP] += e_lj;
253 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Coul] += dvdl_q;
254 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += dvdl_lj;
258 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
260 wallcycle_stop(wcycle, WallCycleCounter::PpPmeWaitRecvF);
263 static void print_large_forces(FILE* fp,
268 ArrayRef<const RVec> x,
269 ArrayRef<const RVec> f)
271 real force2Tolerance = gmx::square(forceTolerance);
272 gmx::index numNonFinite = 0;
273 for (int i = 0; i < md->homenr; i++)
275 real force2 = norm2(f[i]);
276 bool nonFinite = !std::isfinite(force2);
277 if (force2 >= force2Tolerance || nonFinite)
280 "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
293 if (numNonFinite > 0)
295 /* Note that with MPI this fatal call on one rank might interrupt
296 * the printing on other ranks. But we can only avoid that with
297 * an expensive MPI barrier that we would need at each step.
299 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
303 //! When necessary, spreads forces on vsites and computes the virial for \p forceOutputs->forceWithShiftForces()
304 static void postProcessForceWithShiftForces(t_nrnb* nrnb,
305 gmx_wallcycle* wcycle,
307 ArrayRef<const RVec> x,
308 ForceOutputs* forceOutputs,
310 const t_mdatoms& mdatoms,
311 const t_forcerec& fr,
312 gmx::VirtualSitesHandler* vsite,
313 const StepWorkload& stepWork)
315 ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
317 /* If we have NoVirSum forces, but we do not calculate the virial,
318 * we later sum the forceWithShiftForces buffer together with
319 * the noVirSum buffer and spread the combined vsite forces at once.
321 if (vsite && (!forceOutputs->haveForceWithVirial() || stepWork.computeVirial))
323 using VirialHandling = gmx::VirtualSitesHandler::VirialHandling;
325 auto f = forceWithShiftForces.force();
326 auto fshift = forceWithShiftForces.shiftForces();
327 const VirialHandling virialHandling =
328 (stepWork.computeVirial ? VirialHandling::Pbc : VirialHandling::None);
329 vsite->spreadForces(x, f, virialHandling, fshift, nullptr, nrnb, box, wcycle);
330 forceWithShiftForces.haveSpreadVsiteForces() = true;
333 if (stepWork.computeVirial)
335 /* Calculation of the virial must be done after vsites! */
337 0, mdatoms.homenr, as_rvec_array(x.data()), forceWithShiftForces, vir_force, box, nrnb, &fr, fr.pbcType);
341 //! Spread, compute virial for and sum forces, when necessary
342 static void postProcessForces(const t_commrec* cr,
345 gmx_wallcycle* wcycle,
347 ArrayRef<const RVec> x,
348 ForceOutputs* forceOutputs,
350 const t_mdatoms* mdatoms,
351 const t_forcerec* fr,
352 gmx::VirtualSitesHandler* vsite,
353 const StepWorkload& stepWork)
355 // Extract the final output force buffer, which is also the buffer for forces with shift forces
356 ArrayRef<RVec> f = forceOutputs->forceWithShiftForces().force();
358 if (forceOutputs->haveForceWithVirial())
360 auto& forceWithVirial = forceOutputs->forceWithVirial();
364 /* Spread the mesh force on virtual sites to the other particles...
365 * This is parallellized. MPI communication is performed
366 * if the constructing atoms aren't local.
368 GMX_ASSERT(!stepWork.computeVirial || f.data() != forceWithVirial.force_.data(),
369 "We need separate force buffers for shift and virial forces when "
370 "computing the virial");
371 GMX_ASSERT(!stepWork.computeVirial
372 || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
373 "We should spread the force with shift forces separately when computing "
375 const gmx::VirtualSitesHandler::VirialHandling virialHandling =
376 (stepWork.computeVirial ? gmx::VirtualSitesHandler::VirialHandling::NonLinear
377 : gmx::VirtualSitesHandler::VirialHandling::None);
378 matrix virial = { { 0 } };
379 vsite->spreadForces(x, forceWithVirial.force_, virialHandling, {}, virial, nrnb, box, wcycle);
380 forceWithVirial.addVirialContribution(virial);
383 if (stepWork.computeVirial)
385 /* Now add the forces, this is local */
386 sum_forces(f, forceWithVirial.force_);
388 /* Add the direct virial contributions */
390 forceWithVirial.computeVirial_,
391 "forceWithVirial should request virial computation when we request the virial");
392 m_add(vir_force, forceWithVirial.getVirial(), vir_force);
396 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
402 GMX_ASSERT(vsite == nullptr || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
403 "We should have spread the vsite forces (earlier)");
406 if (fr->print_force >= 0)
408 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
412 static void do_nb_verlet(t_forcerec* fr,
413 const interaction_const_t* ic,
414 gmx_enerdata_t* enerd,
415 const StepWorkload& stepWork,
416 const InteractionLocality ilocality,
420 gmx_wallcycle* wcycle)
422 if (!stepWork.computeNonbondedForces)
424 /* skip non-bonded calculation */
428 nonbonded_verlet_t* nbv = fr->nbv.get();
430 /* GPU kernel launch overhead is already timed separately */
433 /* When dynamic pair-list pruning is requested, we need to prune
434 * at nstlistPrune steps.
436 if (nbv->isDynamicPruningStepCpu(step))
438 /* Prune the pair-list beyond fr->ic->rlistPrune using
439 * the current coordinates of the atoms.
441 wallcycle_sub_start(wcycle, WallCycleSubCounter::NonbondedPruning);
442 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
443 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NonbondedPruning);
447 nbv->dispatchNonbondedKernel(
453 enerd->grpp.energyGroupPairTerms[fr->haveBuckingham ? NonBondedEnergyTerms::BuckinghamSR
454 : NonBondedEnergyTerms::LJSR],
455 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR],
459 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
461 int nth = gmx_omp_nthreads_get_simple_rvec_task(ModuleMultiThread::Default, v.ssize());
463 /* Note that we would like to avoid this conditional by putting it
464 * into the omp pragma instead, but then we still take the full
465 * omp parallel for overhead (at least with gcc5).
467 if (!useOpenmpThreading || nth == 1)
476 #pragma omp parallel for num_threads(nth) schedule(static)
477 for (gmx::index i = 0; i < v.ssize(); i++)
484 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
486 * \param groupOptions Group options, containing T-coupling options
488 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
490 real nrdfCoupled = 0;
491 real nrdfUncoupled = 0;
492 real kineticEnergy = 0;
493 for (int g = 0; g < groupOptions.ngtc; g++)
495 if (groupOptions.tau_t[g] >= 0)
497 nrdfCoupled += groupOptions.nrdf[g];
498 kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * gmx::c_boltz;
502 nrdfUncoupled += groupOptions.nrdf[g];
506 /* This conditional with > also catches nrdf=0 */
507 if (nrdfCoupled > nrdfUncoupled)
509 return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
517 /*! \brief This routine checks that the potential energy is finite.
519 * Always checks that the potential energy is finite. If step equals
520 * inputrec.init_step also checks that the magnitude of the potential energy
521 * is reasonable. Terminates with a fatal error when a check fails.
522 * Note that passing this check does not guarantee finite forces,
523 * since those use slightly different arithmetics. But in most cases
524 * there is just a narrow coordinate range where forces are not finite
525 * and energies are finite.
527 * \param[in] step The step number, used for checking and printing
528 * \param[in] enerd The energy data; the non-bonded group energies need to be added to
529 * enerd.term[F_EPOT] before calling this routine \param[in] inputrec The input record
531 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
533 /* Threshold valid for comparing absolute potential energy against
534 * the kinetic energy. Normally one should not consider absolute
535 * potential energy values, but with a factor of one million
536 * we should never get false positives.
538 constexpr real c_thresholdFactor = 1e6;
540 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
541 real averageKineticEnergy = 0;
542 /* We only check for large potential energy at the initial step,
543 * because that is by far the most likely step for this too occur
544 * and because computing the average kinetic energy is not free.
545 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
546 * before they become NaN.
548 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
550 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
553 if (energyIsNotFinite
554 || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
559 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
560 "contributions to the energy are %g and %g, respectively. A %s potential energy "
561 "can be caused by overlapping interactions in bonded interactions or very large%s "
562 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
563 "configuration, incorrect interactions or parameters in the topology.",
566 energyIsNotFinite ? "not finite" : "extremely high",
568 enerd.term[F_COUL_SR],
569 energyIsNotFinite ? "non-finite" : "very high",
570 energyIsNotFinite ? " or Nan" : "");
574 /*! \brief Return true if there are special forces computed this step.
576 * The conditionals exactly correspond to those in computeSpecialForces().
578 static bool haveSpecialForces(const t_inputrec& inputrec,
579 const gmx::ForceProviders& forceProviders,
580 const pull_t* pull_work,
581 const bool computeForces,
585 return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
586 (inputrec.bPull && pull_have_potential(*pull_work)) || // pull
587 inputrec.bRot || // enforced rotation
588 (ed != nullptr) || // flooding
589 (inputrec.bIMD && computeForces)); // IMD
592 /*! \brief Compute forces and/or energies for special algorithms
594 * The intention is to collect all calls to algorithms that compute
595 * forces on local atoms only and that do not contribute to the local
596 * virial sum (but add their virial contribution separately).
597 * Eventually these should likely all become ForceProviders.
598 * Within this function the intention is to have algorithms that do
599 * global communication at the end, so global barriers within the MD loop
600 * are as close together as possible.
602 * \param[in] fplog The log file
603 * \param[in] cr The communication record
604 * \param[in] inputrec The input record
605 * \param[in] awh The Awh module (nullptr if none in use).
606 * \param[in] enforcedRotation Enforced rotation module.
607 * \param[in] imdSession The IMD session
608 * \param[in] pull_work The pull work structure.
609 * \param[in] step The current MD step
610 * \param[in] t The current time
611 * \param[in,out] wcycle Wallcycle accounting struct
612 * \param[in,out] forceProviders Pointer to a list of force providers
613 * \param[in] box The unit cell
614 * \param[in] x The coordinates
615 * \param[in] mdatoms Per atom properties
616 * \param[in] lambda Array of free-energy lambda values
617 * \param[in] stepWork Step schedule flags
618 * \param[in,out] forceWithVirialMtsLevel0 Force and virial for MTS level0 forces
619 * \param[in,out] forceWithVirialMtsLevel1 Force and virial for MTS level1 forces, can be nullptr
620 * \param[in,out] enerd Energy buffer
621 * \param[in,out] ed Essential dynamics pointer
622 * \param[in] didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
624 * \todo Remove didNeighborSearch, which is used incorrectly.
625 * \todo Convert all other algorithms called here to ForceProviders.
627 static void computeSpecialForces(FILE* fplog,
629 const t_inputrec& inputrec,
631 gmx_enfrot* enforcedRotation,
632 gmx::ImdSession* imdSession,
636 gmx_wallcycle* wcycle,
637 gmx::ForceProviders* forceProviders,
639 gmx::ArrayRef<const gmx::RVec> x,
640 const t_mdatoms* mdatoms,
641 gmx::ArrayRef<const real> lambda,
642 const StepWorkload& stepWork,
643 gmx::ForceWithVirial* forceWithVirialMtsLevel0,
644 gmx::ForceWithVirial* forceWithVirialMtsLevel1,
645 gmx_enerdata_t* enerd,
647 bool didNeighborSearch)
649 /* NOTE: Currently all ForceProviders only provide forces.
650 * When they also provide energies, remove this conditional.
652 if (stepWork.computeForces)
654 gmx::ForceProviderInput forceProviderInput(
657 gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->homenr),
658 gmx::arrayRefFromArray(mdatoms->massT, mdatoms->homenr),
662 gmx::ForceProviderOutput forceProviderOutput(forceWithVirialMtsLevel0, enerd);
664 /* Collect forces from modules */
665 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
668 if (inputrec.bPull && pull_have_potential(*pull_work))
670 const int mtsLevel = forceGroupMtsLevel(inputrec.mtsLevels, gmx::MtsForceGroups::Pull);
671 if (mtsLevel == 0 || stepWork.computeSlowForces)
673 auto& forceWithVirial = (mtsLevel == 0) ? forceWithVirialMtsLevel0 : forceWithVirialMtsLevel1;
674 pull_potential_wrapper(
675 cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work, lambda.data(), t, wcycle);
680 const int mtsLevel = forceGroupMtsLevel(inputrec.mtsLevels, gmx::MtsForceGroups::Pull);
681 if (mtsLevel == 0 || stepWork.computeSlowForces)
683 const bool needForeignEnergyDifferences = awh->needForeignEnergyDifferences(step);
684 std::vector<double> foreignLambdaDeltaH, foreignLambdaDhDl;
685 if (needForeignEnergyDifferences)
687 enerd->foreignLambdaTerms.finalizePotentialContributions(
688 enerd->dvdl_lin, lambda, *inputrec.fepvals);
689 std::tie(foreignLambdaDeltaH, foreignLambdaDhDl) = enerd->foreignLambdaTerms.getTerms(cr);
692 auto& forceWithVirial = (mtsLevel == 0) ? forceWithVirialMtsLevel0 : forceWithVirialMtsLevel1;
693 enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(
695 gmx::arrayRefFromArray(mdatoms->massT, mdatoms->nr),
706 /* Add the forces from enforced rotation potentials (if any) */
709 wallcycle_start(wcycle, WallCycleCounter::RotAdd);
710 enerd->term[F_COM_PULL] +=
711 add_rot_forces(enforcedRotation, forceWithVirialMtsLevel0->force_, cr, step, t);
712 wallcycle_stop(wcycle, WallCycleCounter::RotAdd);
717 /* Note that since init_edsam() is called after the initialization
718 * of forcerec, edsam doesn't request the noVirSum force buffer.
719 * Thus if no other algorithm (e.g. PME) requires it, the forces
720 * here will contribute to the virial.
722 do_flood(cr, inputrec, x, forceWithVirialMtsLevel0->force_, ed, box, step, didNeighborSearch);
725 /* Add forces from interactive molecular dynamics (IMD), if any */
726 if (inputrec.bIMD && stepWork.computeForces)
728 imdSession->applyForces(forceWithVirialMtsLevel0->force_);
732 /*! \brief Launch the prepare_step and spread stages of PME GPU.
734 * \param[in] pmedata The PME structure
735 * \param[in] box The box matrix
736 * \param[in] stepWork Step schedule flags
737 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in the device memory.
738 * \param[in] lambdaQ The Coulomb lambda of the current state.
739 * \param[in] wcycle The wallcycle structure
741 static inline void launchPmeGpuSpread(gmx_pme_t* pmedata,
743 const StepWorkload& stepWork,
744 GpuEventSynchronizer* xReadyOnDevice,
746 gmx_wallcycle* wcycle)
748 pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
749 pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle, lambdaQ);
752 /*! \brief Launch the FFT and gather stages of PME GPU
754 * This function only implements setting the output forces (no accumulation).
756 * \param[in] pmedata The PME structure
757 * \param[in] lambdaQ The Coulomb lambda of the current system state.
758 * \param[in] wcycle The wallcycle structure
759 * \param[in] stepWork Step schedule flags
761 static void launchPmeGpuFftAndGather(gmx_pme_t* pmedata,
763 gmx_wallcycle* wcycle,
764 const gmx::StepWorkload& stepWork)
766 pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
767 pme_gpu_launch_gather(pmedata, wcycle, lambdaQ);
771 * Polling wait for either of the PME or nonbonded GPU tasks.
773 * Instead of a static order in waiting for GPU tasks, this function
774 * polls checking which of the two tasks completes first, and does the
775 * associated force buffer reduction overlapped with the other task.
776 * By doing that, unlike static scheduling order, it can always overlap
777 * one of the reductions, regardless of the GPU task completion order.
779 * \param[in] nbv Nonbonded verlet structure
780 * \param[in,out] pmedata PME module data
781 * \param[in,out] forceOutputsNonbonded Force outputs for the non-bonded forces and shift forces
782 * \param[in,out] forceOutputsPme Force outputs for the PME forces and virial
783 * \param[in,out] enerd Energy data structure results are reduced into
784 * \param[in] lambdaQ The Coulomb lambda of the current system state.
785 * \param[in] stepWork Step schedule flags
786 * \param[in] wcycle The wallcycle structure
788 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
790 gmx::ForceOutputs* forceOutputsNonbonded,
791 gmx::ForceOutputs* forceOutputsPme,
792 gmx_enerdata_t* enerd,
794 const StepWorkload& stepWork,
795 gmx_wallcycle* wcycle)
797 bool isPmeGpuDone = false;
798 bool isNbGpuDone = false;
800 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
802 while (!isPmeGpuDone || !isNbGpuDone)
806 GpuTaskCompletion completionType =
807 (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
808 isPmeGpuDone = pme_gpu_try_finish_task(
809 pmedata, stepWork, wcycle, &forceOutputsPme->forceWithVirial(), enerd, lambdaQ, completionType);
814 auto& forceBuffersNonbonded = forceOutputsNonbonded->forceWithShiftForces();
815 GpuTaskCompletion completionType =
816 (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
817 isNbGpuDone = Nbnxm::gpu_try_finish_task(
821 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
822 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
823 forceBuffersNonbonded.shiftForces(),
829 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceBuffersNonbonded.force());
835 /*! \brief Set up the different force buffers; also does clearing.
837 * \param[in] forceHelperBuffers Helper force buffers
838 * \param[in] force force array
839 * \param[in] domainWork Domain lifetime workload flags
840 * \param[in] stepWork Step schedule flags
841 * \param[in] havePpDomainDecomposition Whether we have a PP domain decomposition
842 * \param[out] wcycle wallcycle recording structure
844 * \returns Cleared force output structure
846 static ForceOutputs setupForceOutputs(ForceHelperBuffers* forceHelperBuffers,
847 gmx::ArrayRefWithPadding<gmx::RVec> force,
848 const DomainLifetimeWorkload& domainWork,
849 const StepWorkload& stepWork,
850 const bool havePpDomainDecomposition,
851 gmx_wallcycle* wcycle)
853 wallcycle_sub_start(wcycle, WallCycleSubCounter::ClearForceBuffer);
855 /* NOTE: We assume fr->shiftForces is all zeros here */
856 gmx::ForceWithShiftForces forceWithShiftForces(
857 force, stepWork.computeVirial, forceHelperBuffers->shiftForces());
859 if (stepWork.computeForces
860 && (domainWork.haveCpuLocalForceWork || !stepWork.useGpuFBufferOps
861 || (havePpDomainDecomposition && !stepWork.useGpuFHalo)))
863 /* Clear the short- and long-range forces */
864 clearRVecs(forceWithShiftForces.force(), true);
866 /* Clear the shift forces */
867 clearRVecs(forceWithShiftForces.shiftForces(), false);
870 /* If we need to compute the virial, we might need a separate
871 * force buffer for algorithms for which the virial is calculated
872 * directly, such as PME. Otherwise, forceWithVirial uses the
873 * the same force (f in legacy calls) buffer as other algorithms.
875 const bool useSeparateForceWithVirialBuffer =
876 (stepWork.computeForces
877 && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
878 /* forceWithVirial uses the local atom range only */
879 gmx::ForceWithVirial forceWithVirial(
880 useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
881 : force.unpaddedArrayRef(),
882 stepWork.computeVirial);
884 if (useSeparateForceWithVirialBuffer)
886 /* TODO: update comment
887 * We only compute forces on local atoms. Note that vsites can
888 * spread to non-local atoms, but that part of the buffer is
889 * cleared separately in the vsite spreading code.
891 clearRVecs(forceWithVirial.force_, true);
894 wallcycle_sub_stop(wcycle, WallCycleSubCounter::ClearForceBuffer);
897 forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(), forceWithVirial);
901 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
903 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec& inputrec,
904 const t_forcerec& fr,
905 const pull_t* pull_work,
907 const t_mdatoms& mdatoms,
908 const SimulationWorkload& simulationWork,
909 const StepWorkload& stepWork)
911 DomainLifetimeWorkload domainWork;
912 // Note that haveSpecialForces is constant over the whole run
913 domainWork.haveSpecialForces =
914 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
915 domainWork.haveCpuListedForceWork = false;
916 domainWork.haveCpuBondedWork = false;
917 for (const auto& listedForces : fr.listedForces)
919 if (listedForces.haveCpuListedForces(*fr.fcdata))
921 domainWork.haveCpuListedForceWork = true;
923 if (listedForces.haveCpuBondeds())
925 domainWork.haveCpuBondedWork = true;
928 domainWork.haveGpuBondedWork =
929 ((fr.listedForcesGpu != nullptr) && fr.listedForcesGpu->haveInteractions());
930 // Note that haveFreeEnergyWork is constant over the whole run
931 domainWork.haveFreeEnergyWork =
932 (fr.efep != FreeEnergyPerturbationType::No && mdatoms.nPerturbed != 0);
933 // We assume we have local force work if there are CPU
934 // force tasks including PME or nonbondeds.
935 domainWork.haveCpuLocalForceWork =
936 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
937 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
938 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
939 domainWork.haveLocalForceContribInCpuBuffer =
940 domainWork.haveCpuLocalForceWork || simulationWork.havePpDomainDecomposition;
941 domainWork.haveNonLocalForceContribInCpuBuffer =
942 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
947 /*! \brief Set up force flag stuct from the force bitmask.
949 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
950 * \param[in] mtsLevels The multiple time-stepping levels, either empty or 2 levels
951 * \param[in] step The current MD step
952 * \param[in] simulationWork Simulation workload description.
953 * \param[in] rankHasPmeDuty If this rank computes PME.
955 * \returns New Stepworkload description.
957 static StepWorkload setupStepWorkload(const int legacyFlags,
958 ArrayRef<const gmx::MtsLevel> mtsLevels,
960 const SimulationWorkload& simulationWork,
961 const bool rankHasPmeDuty)
963 GMX_ASSERT(mtsLevels.empty() || mtsLevels.size() == 2, "Expect 0 or 2 MTS levels");
964 const bool computeSlowForces = (mtsLevels.empty() || step % mtsLevels[1].stepFactor == 0);
967 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
968 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
969 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
970 flags.computeSlowForces = computeSlowForces;
971 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
972 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
973 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
974 flags.useOnlyMtsCombinedForceBuffer = ((legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0);
975 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
976 flags.computeNonbondedForces =
977 ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded
978 && !(simulationWork.computeNonbondedAtMtsLevel1 && !computeSlowForces);
979 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
981 if (simulationWork.useGpuBufferOps)
983 GMX_ASSERT(simulationWork.useGpuNonbonded,
984 "Can only offload buffer ops if nonbonded computation is also offloaded");
986 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
987 // on virial steps the CPU reduction path is taken
988 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
989 flags.useGpuPmeFReduction = flags.computeSlowForces && flags.useGpuFBufferOps && simulationWork.useGpuPme
990 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication);
991 flags.useGpuXHalo = simulationWork.useGpuHaloExchange;
992 flags.useGpuFHalo = simulationWork.useGpuHaloExchange && flags.useGpuFBufferOps;
993 // Note that rankHasPmeDuty is used confusingly due to the way cr->duty is set up (can be true even for non-PME runs),
994 // but the haveGpuPmeOnThisRank still ends up correct as simulationWork.useGpuPme == false in such cases.
995 // TODO: improve this when duty-reliance is eliminated
996 flags.haveGpuPmeOnThisRank = simulationWork.useGpuPme && rankHasPmeDuty && flags.computeSlowForces;
997 flags.combineMtsForcesBeforeHaloExchange =
998 (flags.computeForces && simulationWork.useMts && flags.computeSlowForces
999 && flags.useOnlyMtsCombinedForceBuffer
1000 && !(flags.computeVirial || simulationWork.useGpuNonbonded || flags.haveGpuPmeOnThisRank));
1006 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
1009 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
1010 gmx::ListedForcesGpu* listedForcesGpu,
1012 gmx_enerdata_t* enerd,
1013 const gmx::MdrunScheduleWorkload& runScheduleWork,
1015 gmx_wallcycle* wcycle)
1017 if (runScheduleWork.simulationWork.useGpuNonbonded && runScheduleWork.stepWork.computeNonbondedForces)
1019 /* Launch pruning before buffer clearing because the API overhead of the
1020 * clear kernel launches can leave the GPU idle while it could be running
1023 if (nbv->isDynamicPruningStepGpu(step))
1025 nbv->dispatchPruneKernelGpu(step);
1028 /* now clear the GPU outputs while we finish the step on the CPU */
1029 wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1030 wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1031 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
1032 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1033 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1036 if (runScheduleWork.stepWork.haveGpuPmeOnThisRank)
1038 pme_gpu_reinit_computation(pmedata, wcycle);
1041 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
1043 // in principle this should be included in the DD balancing region,
1044 // but generally it is infrequent so we'll omit it for the sake of
1046 listedForcesGpu->waitAccumulateEnergyTerms(enerd);
1048 listedForcesGpu->clearEnergies();
1052 //! \brief Data structure to hold dipole-related data and staging arrays
1055 //! Dipole staging for fast summing over MPI
1056 gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
1057 //! Dipole staging for states A and B (index 0 and 1 resp.)
1058 gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
1062 static void reduceAndUpdateMuTot(DipoleData* dipoleData,
1063 const t_commrec* cr,
1064 const bool haveFreeEnergy,
1065 gmx::ArrayRef<const real> lambda,
1067 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1071 gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
1072 ddBalanceRegionHandler.reopenRegionCpu();
1074 for (int i = 0; i < 2; i++)
1076 for (int j = 0; j < DIM; j++)
1078 dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
1082 if (!haveFreeEnergy)
1084 copy_rvec(dipoleData->muStateAB[0], muTotal);
1088 for (int j = 0; j < DIM; j++)
1090 muTotal[j] = (1.0 - lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)])
1091 * dipoleData->muStateAB[0][j]
1092 + lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]
1093 * dipoleData->muStateAB[1][j];
1098 /*! \brief Combines MTS level0 and level1 force buffers into a full and MTS-combined force buffer.
1100 * \param[in] numAtoms The number of atoms to combine forces for
1101 * \param[in,out] forceMtsLevel0 Input: F_level0, output: F_level0 + F_level1
1102 * \param[in,out] forceMts Input: F_level1, output: F_level0 + mtsFactor * F_level1
1103 * \param[in] mtsFactor The factor between the level0 and level1 time step
1105 static void combineMtsForces(const int numAtoms,
1106 ArrayRef<RVec> forceMtsLevel0,
1107 ArrayRef<RVec> forceMts,
1108 const real mtsFactor)
1110 const int gmx_unused numThreads = gmx_omp_nthreads_get(ModuleMultiThread::Default);
1111 #pragma omp parallel for num_threads(numThreads) schedule(static)
1112 for (int i = 0; i < numAtoms; i++)
1114 const RVec forceMtsLevel0Tmp = forceMtsLevel0[i];
1115 forceMtsLevel0[i] += forceMts[i];
1116 forceMts[i] = forceMtsLevel0Tmp + mtsFactor * forceMts[i];
1120 /*! \brief Setup for the local and non-local GPU force reductions:
1121 * reinitialization plus the registration of forces and dependencies.
1123 * \param [in] runScheduleWork Schedule workload flag structure
1124 * \param [in] cr Communication record object
1125 * \param [in] fr Force record object
1127 static void setupGpuForceReductions(gmx::MdrunScheduleWorkload* runScheduleWork,
1128 const t_commrec* cr,
1132 nonbonded_verlet_t* nbv = fr->nbv.get();
1133 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1135 // (re-)initialize local GPU force reduction
1136 const bool accumulate = runScheduleWork->domainWork.haveCpuLocalForceWork
1137 || runScheduleWork->simulationWork.havePpDomainDecomposition;
1138 const int atomStart = 0;
1139 fr->gpuForceReduction[gmx::AtomLocality::Local]->reinit(stateGpu->getForces(),
1140 nbv->getNumAtoms(AtomLocality::Local),
1141 nbv->getGridIndices(),
1144 stateGpu->fReducedOnDevice());
1146 // register forces and add dependencies
1147 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerNbnxmForce(Nbnxm::gpu_get_f(nbv->gpu_nbv));
1149 if (runScheduleWork->simulationWork.useGpuPme
1150 && (!runScheduleWork->simulationWork.haveSeparatePmeRank
1151 || runScheduleWork->simulationWork.useGpuPmePpCommunication))
1153 DeviceBuffer<gmx::RVec> forcePtr =
1154 runScheduleWork->simulationWork.haveSeparatePmeRank
1155 ? fr->pmePpCommGpu->getGpuForceStagingPtr() // buffer received from other GPU
1156 : pme_gpu_get_device_f(fr->pmedata); // PME force buffer on same GPU
1157 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerRvecForce(forcePtr);
1159 GpuEventSynchronizer* const pmeSynchronizer =
1160 (runScheduleWork->simulationWork.haveSeparatePmeRank
1161 ? fr->pmePpCommGpu->getForcesReadySynchronizer() // buffer received from other GPU
1162 : pme_gpu_get_f_ready_synchronizer(fr->pmedata)); // PME force buffer on same GPU
1165 GMX_ASSERT(pmeSynchronizer != nullptr, "PME force ready cuda event should not be NULL");
1166 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(pmeSynchronizer);
1170 if (runScheduleWork->domainWork.haveCpuLocalForceWork && !runScheduleWork->simulationWork.useGpuHaloExchange)
1172 // in the DD case we use the same stream for H2D and reduction, hence no explicit dependency needed
1173 if (!runScheduleWork->simulationWork.havePpDomainDecomposition)
1175 const bool useGpuForceBufferOps = true;
1176 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1177 stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::All, useGpuForceBufferOps));
1181 if (runScheduleWork->simulationWork.useGpuHaloExchange)
1183 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1184 cr->dd->gpuHaloExchange[0][0]->getForcesReadyOnDeviceEvent());
1187 if (runScheduleWork->simulationWork.havePpDomainDecomposition)
1189 // (re-)initialize non-local GPU force reduction
1190 const bool accumulate = runScheduleWork->domainWork.haveCpuBondedWork
1191 || runScheduleWork->domainWork.haveFreeEnergyWork;
1192 const int atomStart = dd_numHomeAtoms(*cr->dd);
1193 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->reinit(stateGpu->getForces(),
1194 nbv->getNumAtoms(AtomLocality::NonLocal),
1195 nbv->getGridIndices(),
1199 // register forces and add dependencies
1200 // in the DD case we use the same stream for H2D and reduction, hence no explicit dependency needed
1201 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->registerNbnxmForce(
1202 Nbnxm::gpu_get_f(nbv->gpu_nbv));
1207 /*! \brief Return the number of local atoms.
1209 static int getLocalAtomCount(const gmx_domdec_t* dd, const t_mdatoms& mdatoms, bool havePPDomainDecomposition)
1211 GMX_ASSERT(!(havePPDomainDecomposition && (dd == nullptr)),
1212 "Can't have PP decomposition with dd uninitialized!");
1213 return havePPDomainDecomposition ? dd_numAtomsZones(*dd) : mdatoms.homenr;
1217 void do_force(FILE* fplog,
1218 const t_commrec* cr,
1219 const gmx_multisim_t* ms,
1220 const t_inputrec& inputrec,
1222 gmx_enfrot* enforcedRotation,
1223 gmx::ImdSession* imdSession,
1227 gmx_wallcycle* wcycle,
1228 const gmx_localtop_t* top,
1230 gmx::ArrayRefWithPadding<gmx::RVec> x,
1231 const history_t* hist,
1232 gmx::ForceBuffersView* forceView,
1234 const t_mdatoms* mdatoms,
1235 gmx_enerdata_t* enerd,
1236 gmx::ArrayRef<const real> lambda,
1238 gmx::MdrunScheduleWorkload* runScheduleWork,
1239 gmx::VirtualSitesHandler* vsite,
1244 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1246 auto force = forceView->forceWithPadding();
1247 GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1248 "The size of the force buffer should be at least the number of atoms to compute "
1251 nonbonded_verlet_t* nbv = fr->nbv.get();
1252 interaction_const_t* ic = fr->ic.get();
1254 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1256 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1258 runScheduleWork->stepWork = setupStepWorkload(
1259 legacyFlags, inputrec.mtsLevels, step, simulationWork, thisRankHasDuty(cr, DUTY_PME));
1260 const StepWorkload& stepWork = runScheduleWork->stepWork;
1262 /* At a search step we need to start the first balancing region
1263 * somewhere early inside the step after communication during domain
1264 * decomposition (and not during the previous step as usual).
1266 if (stepWork.doNeighborSearch)
1268 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1271 clear_mat(vir_force);
1273 if (fr->pbcType != PbcType::No)
1275 /* Compute shift vectors every step,
1276 * because of pressure coupling or box deformation!
1278 if (stepWork.haveDynamicBox && stepWork.stateChanged)
1280 calc_shifts(box, fr->shift_vec);
1283 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1284 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1287 put_atoms_in_box_omp(fr->pbcType,
1289 x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1290 gmx_omp_nthreads_get(ModuleMultiThread::Default));
1291 inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1295 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1297 const bool pmeSendCoordinatesFromGpu =
1298 simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1299 const bool reinitGpuPmePpComms =
1300 simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1302 auto* localXReadyOnDevice = (stepWork.haveGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1303 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1304 AtomLocality::Local, simulationWork, stepWork)
1307 GMX_ASSERT(simulationWork.useGpuHaloExchange
1308 == ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange[0].empty())),
1309 "The GPU halo exchange is active, but it has not been constructed.");
1311 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1312 // Copy coordinate from the GPU if update is on the GPU and there
1313 // are forces to be computed on the CPU, or for the computation of
1314 // virial, or if host-side data will be transferred from this task
1315 // to a remote task for halo exchange or PME-PP communication. At
1316 // search steps the current coordinates are already on the host,
1317 // hence copy is not needed.
1318 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1319 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1320 || simulationWork.useCpuPmePpCommunication || simulationWork.useCpuHaloExchange
1321 || simulationWork.computeMuTot))
1323 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1324 haveCopiedXFromGpu = true;
1327 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1328 // The local coordinates can be copied right away.
1329 // NOTE: Consider moving this copy to right after they are updated and constrained,
1330 // if the later is not offloaded.
1331 if (stepWork.haveGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1333 if (stepWork.doNeighborSearch)
1335 // TODO refactor this to do_md, after partitioning.
1336 stateGpu->reinit(mdatoms->homenr,
1337 getLocalAtomCount(cr->dd, *mdatoms, simulationWork.havePpDomainDecomposition));
1338 if (stepWork.haveGpuPmeOnThisRank)
1340 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1341 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1344 // We need to copy coordinates when:
1345 // 1. Update is not offloaded
1346 // 2. The buffers were reinitialized on search step
1347 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1349 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1350 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1354 if (simulationWork.haveSeparatePmeRank && stepWork.computeSlowForces)
1356 /* Send particle coordinates to the pme nodes */
1357 if (!pmeSendCoordinatesFromGpu && !stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1359 GMX_ASSERT(haveCopiedXFromGpu,
1360 "a wait should only be triggered if copy has been scheduled");
1361 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1364 gmx_pme_send_coordinates(fr,
1367 x.unpaddedArrayRef(),
1368 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1369 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1370 (stepWork.computeVirial || stepWork.computeEnergy),
1372 simulationWork.useGpuPmePpCommunication,
1373 reinitGpuPmePpComms,
1374 pmeSendCoordinatesFromGpu,
1375 stepWork.useGpuPmeFReduction,
1376 localXReadyOnDevice,
1380 if (stepWork.haveGpuPmeOnThisRank)
1382 launchPmeGpuSpread(fr->pmedata,
1385 localXReadyOnDevice,
1386 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1390 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1392 /* do gridding for pair search */
1393 if (stepWork.doNeighborSearch)
1395 if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1397 fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1400 wallcycle_start(wcycle, WallCycleCounter::NS);
1401 if (!DOMAINDECOMP(cr))
1403 const rvec vzero = { 0.0_real, 0.0_real, 0.0_real };
1404 const rvec boxDiagonal = { box[XX][XX], box[YY][YY], box[ZZ][ZZ] };
1405 wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSGridLocal);
1406 nbnxn_put_on_grid(nbv,
1412 { 0, mdatoms->homenr },
1415 x.unpaddedArrayRef(),
1418 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSGridLocal);
1422 wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSGridNonLocal);
1423 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->atomInfo, x.unpaddedArrayRef());
1424 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSGridNonLocal);
1427 nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1428 gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1431 wallcycle_stop(wcycle, WallCycleCounter::NS);
1433 /* initialize the GPU nbnxm atom data and bonded data structures */
1434 if (simulationWork.useGpuNonbonded)
1436 // Note: cycle counting only nononbondeds, GPU listed forces counts internally
1437 wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1438 wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1439 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1440 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1441 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1443 if (fr->listedForcesGpu)
1445 /* Now we put all atoms on the grid, we can assign bonded
1446 * interactions to the GPU, where the grid order is
1447 * needed. Also the xq, f and fshift device buffers have
1448 * been reallocated if needed, so the bonded code can
1449 * learn about them. */
1450 // TODO the xq, f, and fshift buffers are now shared
1451 // resources, so they should be maintained by a
1452 // higher-level object than the nb module.
1453 fr->listedForcesGpu->updateInteractionListsAndDeviceBuffers(
1454 nbv->getGridIndices(),
1456 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1457 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1458 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1462 // Need to run after the GPU-offload bonded interaction lists
1463 // are set up to be able to determine whether there is bonded work.
1464 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1465 inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1467 wallcycle_start_nocount(wcycle, WallCycleCounter::NS);
1468 wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSSearchLocal);
1469 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1470 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1472 nbv->setupGpuShortRangeWork(fr->listedForcesGpu.get(), InteractionLocality::Local);
1474 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSSearchLocal);
1475 wallcycle_stop(wcycle, WallCycleCounter::NS);
1477 if (stepWork.useGpuXBufferOps)
1479 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1482 if (simulationWork.useGpuBufferOps)
1484 setupGpuForceReductions(runScheduleWork, cr, fr);
1487 else if (!EI_TPI(inputrec.eI) && stepWork.computeNonbondedForces)
1489 if (stepWork.useGpuXBufferOps)
1491 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1492 nbv->convertCoordinatesGpu(AtomLocality::Local, stateGpu->getCoordinates(), localXReadyOnDevice);
1496 if (simulationWork.useGpuUpdate)
1498 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1499 GMX_ASSERT(haveCopiedXFromGpu,
1500 "a wait should only be triggered if copy has been scheduled");
1501 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1503 nbv->convertCoordinates(AtomLocality::Local, x.unpaddedArrayRef());
1507 if (simulationWork.useGpuNonbonded && (stepWork.computeNonbondedForces || domainWork.haveGpuBondedWork))
1509 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1511 wallcycle_start(wcycle, WallCycleCounter::LaunchGpu);
1512 wallcycle_sub_start(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1513 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1514 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1516 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1518 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1519 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1520 // with X buffer ops offloaded to the GPU on all but the search steps
1522 // bonded work not split into separate local and non-local, so with DD
1523 // we can only launch the kernel after non-local coordinates have been received.
1524 if (domainWork.haveGpuBondedWork && !simulationWork.havePpDomainDecomposition)
1526 fr->listedForcesGpu->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1529 /* launch local nonbonded work on GPU */
1530 wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1531 wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1532 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1533 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1534 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1537 if (stepWork.haveGpuPmeOnThisRank)
1539 // In PME GPU and mixed mode we launch FFT / gather after the
1540 // X copy/transform to allow overlap as well as after the GPU NB
1541 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1542 // the nonbonded kernel.
1543 launchPmeGpuFftAndGather(fr->pmedata,
1544 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1549 /* Communicate coordinates and sum dipole if necessary +
1550 do non-local pair search */
1551 if (simulationWork.havePpDomainDecomposition)
1553 if (stepWork.doNeighborSearch)
1555 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1556 wallcycle_start_nocount(wcycle, WallCycleCounter::NS);
1557 wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSSearchNonLocal);
1558 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1559 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1561 nbv->setupGpuShortRangeWork(fr->listedForcesGpu.get(), InteractionLocality::NonLocal);
1562 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSSearchNonLocal);
1563 wallcycle_stop(wcycle, WallCycleCounter::NS);
1564 // TODO refactor this GPU halo exchange re-initialisation
1565 // to location in do_md where GPU halo exchange is
1566 // constructed at partitioning, after above stateGpu
1567 // re-initialization has similarly been refactored
1568 if (simulationWork.useGpuHaloExchange)
1570 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1575 if (stepWork.useGpuXHalo)
1577 // The following must be called after local setCoordinates (which records an event
1578 // when the coordinate data has been copied to the device).
1579 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1581 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1583 // non-local part of coordinate buffer must be copied back to host for CPU work
1584 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1589 if (simulationWork.useGpuUpdate)
1591 GMX_ASSERT(haveCopiedXFromGpu,
1592 "a wait should only be triggered if copy has been scheduled");
1593 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1595 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1598 if (stepWork.useGpuXBufferOps)
1600 if (!stepWork.haveGpuPmeOnThisRank && !stepWork.useGpuXHalo)
1602 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1604 nbv->convertCoordinatesGpu(AtomLocality::NonLocal,
1605 stateGpu->getCoordinates(),
1606 stateGpu->getCoordinatesReadyOnDeviceEvent(
1607 AtomLocality::NonLocal, simulationWork, stepWork));
1611 nbv->convertCoordinates(AtomLocality::NonLocal, x.unpaddedArrayRef());
1615 if (simulationWork.useGpuNonbonded)
1618 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1620 wallcycle_start(wcycle, WallCycleCounter::LaunchGpu);
1621 wallcycle_sub_start(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1622 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1623 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1624 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1627 if (domainWork.haveGpuBondedWork)
1629 fr->listedForcesGpu->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1632 /* launch non-local nonbonded tasks on GPU */
1633 wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1634 wallcycle_sub_start(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1635 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1636 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1637 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1641 if (simulationWork.useGpuNonbonded && stepWork.computeNonbondedForces)
1643 /* launch D2H copy-back F */
1644 wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1645 wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1647 if (simulationWork.havePpDomainDecomposition)
1649 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1651 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1652 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1654 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1656 fr->listedForcesGpu->launchEnergyTransfer();
1658 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1661 gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1662 if (fr->wholeMoleculeTransform)
1664 xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1667 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1668 // this wait ensures that the D2H transfer is complete.
1669 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch)
1671 const bool needCoordsOnHost = (runScheduleWork->domainWork.haveCpuLocalForceWork
1672 || stepWork.computeVirial || simulationWork.computeMuTot);
1673 const bool haveAlreadyWaited = simulationWork.useCpuHaloExchange;
1674 if (needCoordsOnHost && !haveAlreadyWaited)
1676 GMX_ASSERT(haveCopiedXFromGpu,
1677 "a wait should only be triggered if copy has been scheduled");
1678 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1682 DipoleData dipoleData;
1684 if (simulationWork.computeMuTot)
1686 const int start = 0;
1688 /* Calculate total (local) dipole moment in a temporary common array.
1689 * This makes it possible to sum them over nodes faster.
1691 gmx::ArrayRef<const gmx::RVec> xRef =
1692 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1696 mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
1697 : gmx::ArrayRef<real>{},
1698 mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
1699 : gmx::ArrayRef<real>{},
1700 mdatoms->nChargePerturbed != 0,
1701 dipoleData.muStaging[0],
1702 dipoleData.muStaging[1]);
1704 reduceAndUpdateMuTot(
1705 &dipoleData, cr, (fr->efep != FreeEnergyPerturbationType::No), lambda, muTotal, ddBalanceRegionHandler);
1708 /* Reset energies */
1709 reset_enerdata(enerd);
1711 if (DOMAINDECOMP(cr) && simulationWork.haveSeparatePmeRank)
1713 wallcycle_start(wcycle, WallCycleCounter::PpDuringPme);
1714 dd_force_flop_start(cr->dd, nrnb);
1719 wallcycle_start(wcycle, WallCycleCounter::Rot);
1720 do_rotation(cr, enforcedRotation, box, x.unpaddedConstArrayRef(), t, step, stepWork.doNeighborSearch);
1721 wallcycle_stop(wcycle, WallCycleCounter::Rot);
1724 /* Start the force cycle counter.
1725 * Note that a different counter is used for dynamic load balancing.
1727 wallcycle_start(wcycle, WallCycleCounter::Force);
1729 /* Set up and clear force outputs:
1730 * forceOutMtsLevel0: everything except what is in the other two outputs
1731 * forceOutMtsLevel1: PME-mesh and listed-forces group 1
1732 * forceOutNonbonded: non-bonded forces
1733 * Without multiple time stepping all point to the same object.
1734 * With multiple time-stepping the use is different for MTS fast (level0 only) and slow steps.
1736 ForceOutputs forceOutMtsLevel0 = setupForceOutputs(
1737 &fr->forceHelperBuffers[0], force, domainWork, stepWork, simulationWork.havePpDomainDecomposition, wcycle);
1739 // Force output for MTS combined forces, only set at level1 MTS steps
1740 std::optional<ForceOutputs> forceOutMts =
1741 (simulationWork.useMts && stepWork.computeSlowForces)
1742 ? std::optional(setupForceOutputs(&fr->forceHelperBuffers[1],
1743 forceView->forceMtsCombinedWithPadding(),
1746 simulationWork.havePpDomainDecomposition,
1750 ForceOutputs* forceOutMtsLevel1 =
1751 simulationWork.useMts ? (stepWork.computeSlowForces ? &forceOutMts.value() : nullptr)
1752 : &forceOutMtsLevel0;
1754 const bool nonbondedAtMtsLevel1 = runScheduleWork->simulationWork.computeNonbondedAtMtsLevel1;
1756 ForceOutputs* forceOutNonbonded = nonbondedAtMtsLevel1 ? forceOutMtsLevel1 : &forceOutMtsLevel0;
1758 if (inputrec.bPull && pull_have_constraint(*pull_work))
1760 clear_pull_forces(pull_work);
1763 /* We calculate the non-bonded forces, when done on the CPU, here.
1764 * We do this before calling do_force_lowlevel, because in that
1765 * function, the listed forces are calculated before PME, which
1766 * does communication. With this order, non-bonded and listed
1767 * force calculation imbalance can be balanced out by the domain
1768 * decomposition load balancing.
1771 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1773 if (!useOrEmulateGpuNb)
1775 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1778 if (fr->efep != FreeEnergyPerturbationType::No && stepWork.computeNonbondedForces)
1780 /* Calculate the local and non-local free energy interactions here.
1781 * Happens here on the CPU both with and without GPU.
1783 nbv->dispatchFreeEnergyKernel(
1784 InteractionLocality::Local,
1785 x.unpaddedArrayRef(),
1786 &forceOutNonbonded->forceWithShiftForces(),
1787 fr->use_simd_kernels,
1794 mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
1795 : gmx::ArrayRef<real>{},
1796 mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
1797 : gmx::ArrayRef<real>{},
1798 mdatoms->typeA ? gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr)
1799 : gmx::ArrayRef<int>{},
1800 mdatoms->typeB ? gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr)
1801 : gmx::ArrayRef<int>{},
1802 inputrec.fepvals.get(),
1808 if (simulationWork.havePpDomainDecomposition)
1810 nbv->dispatchFreeEnergyKernel(
1811 InteractionLocality::NonLocal,
1812 x.unpaddedArrayRef(),
1813 &forceOutNonbonded->forceWithShiftForces(),
1814 fr->use_simd_kernels,
1821 mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
1822 : gmx::ArrayRef<real>{},
1823 mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
1824 : gmx::ArrayRef<real>{},
1825 mdatoms->typeA ? gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr)
1826 : gmx::ArrayRef<int>{},
1827 mdatoms->typeB ? gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr)
1828 : gmx::ArrayRef<int>{},
1829 inputrec.fepvals.get(),
1837 if (stepWork.computeNonbondedForces && !useOrEmulateGpuNb)
1839 if (simulationWork.havePpDomainDecomposition)
1841 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1844 if (stepWork.computeForces)
1846 /* Add all the non-bonded force to the normal force array.
1847 * This can be split into a local and a non-local part when overlapping
1848 * communication with calculation with domain decomposition.
1850 wallcycle_stop(wcycle, WallCycleCounter::Force);
1851 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All,
1852 forceOutNonbonded->forceWithShiftForces().force());
1853 wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
1856 /* If there are multiple fshift output buffers we need to reduce them */
1857 if (stepWork.computeVirial)
1859 /* This is not in a subcounter because it takes a
1860 negligible and constant-sized amount of time */
1861 nbnxn_atomdata_add_nbat_fshift_to_fshift(
1862 *nbv->nbat, forceOutNonbonded->forceWithShiftForces().shiftForces());
1866 // TODO Force flags should include haveFreeEnergyWork for this domain
1867 if (stepWork.useGpuXHalo && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1869 wallcycle_stop(wcycle, WallCycleCounter::Force);
1870 /* Wait for non-local coordinate data to be copied from device */
1871 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1872 wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
1875 // Compute wall interactions, when present.
1876 // Note: should be moved to special forces.
1877 if (inputrec.nwall && stepWork.computeNonbondedForces)
1879 /* foreign lambda component for walls */
1880 real dvdl_walls = do_walls(inputrec,
1883 mdatoms->typeA ? gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr)
1884 : gmx::ArrayRef<int>{},
1885 mdatoms->typeB ? gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr)
1886 : gmx::ArrayRef<int>{},
1887 mdatoms->cENER ? gmx::arrayRefFromArray(mdatoms->cENER, mdatoms->nr)
1888 : gmx::ArrayRef<unsigned short>{},
1890 mdatoms->nPerturbed,
1891 x.unpaddedConstArrayRef(),
1892 &forceOutMtsLevel0.forceWithVirial(),
1893 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1894 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR],
1896 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += dvdl_walls;
1899 if (stepWork.computeListedForces)
1901 /* Check whether we need to take into account PBC in listed interactions */
1902 bool needMolPbc = false;
1903 for (const auto& listedForces : fr->listedForces)
1905 if (listedForces.haveCpuListedForces(*fr->fcdata))
1907 needMolPbc = fr->bMolPBC;
1915 /* Since all atoms are in the rectangular or triclinic unit-cell,
1916 * only single box vector shifts (2 in x) are required.
1918 set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1921 for (int mtsIndex = 0; mtsIndex < (simulationWork.useMts && stepWork.computeSlowForces ? 2 : 1);
1924 ListedForces& listedForces = fr->listedForces[mtsIndex];
1925 ForceOutputs& forceOut = (mtsIndex == 0 ? forceOutMtsLevel0 : *forceOutMtsLevel1);
1926 listedForces.calculate(wcycle,
1928 inputrec.fepvals.get(),
1942 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr,
1947 if (stepWork.computeSlowForces)
1949 calculateLongRangeNonbondeds(fr,
1955 x.unpaddedConstArrayRef(),
1956 &forceOutMtsLevel1->forceWithVirial(),
1960 dipoleData.muStateAB,
1962 ddBalanceRegionHandler);
1965 wallcycle_stop(wcycle, WallCycleCounter::Force);
1967 // VdW dispersion correction, only computed on master rank to avoid double counting
1968 if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1970 // Calculate long range corrections to pressure and energy
1971 const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
1972 box, lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)]);
1974 if (stepWork.computeEnergy)
1976 enerd->term[F_DISPCORR] = correction.energy;
1977 enerd->term[F_DVDL_VDW] += correction.dvdl;
1978 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += correction.dvdl;
1980 if (stepWork.computeVirial)
1982 correction.correctVirial(vir_force);
1983 enerd->term[F_PDISPCORR] = correction.pressure;
1987 computeSpecialForces(fplog,
1999 x.unpaddedArrayRef(),
2003 &forceOutMtsLevel0.forceWithVirial(),
2004 forceOutMtsLevel1 ? &forceOutMtsLevel1->forceWithVirial() : nullptr,
2007 stepWork.doNeighborSearch);
2009 if (simulationWork.havePpDomainDecomposition && stepWork.computeForces && stepWork.useGpuFHalo
2010 && domainWork.haveCpuLocalForceWork)
2012 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(), AtomLocality::Local);
2015 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
2016 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
2017 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFHalo),
2018 "The schedule below does not allow for nonbonded MTS with GPU halo exchange");
2019 // Will store the amount of cycles spent waiting for the GPU that
2020 // will be later used in the DLB accounting.
2021 float cycles_wait_gpu = 0;
2022 if (useOrEmulateGpuNb && stepWork.computeNonbondedForces)
2024 auto& forceWithShiftForces = forceOutNonbonded->forceWithShiftForces();
2026 /* wait for non-local forces (or calculate in emulation mode) */
2027 if (simulationWork.havePpDomainDecomposition)
2029 if (simulationWork.useGpuNonbonded)
2031 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
2034 AtomLocality::NonLocal,
2035 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
2036 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
2037 forceWithShiftForces.shiftForces(),
2042 wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
2044 fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes, step, nrnb, wcycle);
2045 wallcycle_stop(wcycle, WallCycleCounter::Force);
2048 if (stepWork.useGpuFBufferOps)
2050 if (domainWork.haveNonLocalForceContribInCpuBuffer)
2052 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2053 AtomLocality::NonLocal);
2057 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->execute();
2059 if (!stepWork.useGpuFHalo)
2061 // copy from GPU input for dd_move_f()
2062 stateGpu->copyForcesFromGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2063 AtomLocality::NonLocal);
2068 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
2071 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
2073 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
2078 /* Combining the forces for multiple time stepping before the halo exchange, when possible,
2079 * avoids an extra halo exchange (when DD is used) and post-processing step.
2081 if (stepWork.combineMtsForcesBeforeHaloExchange)
2083 combineMtsForces(getLocalAtomCount(cr->dd, *mdatoms, simulationWork.havePpDomainDecomposition),
2084 force.unpaddedArrayRef(),
2085 forceView->forceMtsCombined(),
2086 inputrec.mtsLevels[1].stepFactor);
2089 if (simulationWork.havePpDomainDecomposition)
2091 /* We are done with the CPU compute.
2092 * We will now communicate the non-local forces.
2093 * If we use a GPU this will overlap with GPU work, so in that case
2094 * we do not close the DD force balancing region here.
2096 ddBalanceRegionHandler.closeAfterForceComputationCpu();
2098 if (stepWork.computeForces)
2101 if (stepWork.useGpuFHalo)
2103 // If there exist CPU forces, data from halo exchange should accumulate into these
2104 bool accumulateForces = domainWork.haveCpuLocalForceWork;
2105 if (!accumulateForces)
2107 // Force halo exchange will set a subset of local atoms with remote non-local data
2108 // First clear local portion of force array, so that untouched atoms are zero
2109 stateGpu->clearForcesOnGpu(AtomLocality::Local);
2111 communicateGpuHaloForces(*cr, accumulateForces);
2115 if (stepWork.useGpuFBufferOps)
2117 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
2120 // Without MTS or with MTS at slow steps with uncombined forces we need to
2121 // communicate the fast forces
2122 if (!simulationWork.useMts || !stepWork.combineMtsForcesBeforeHaloExchange)
2124 dd_move_f(cr->dd, &forceOutMtsLevel0.forceWithShiftForces(), wcycle);
2126 // With MTS we need to communicate the slow or combined (in forceOutMtsLevel1) forces
2127 if (simulationWork.useMts && stepWork.computeSlowForces)
2129 dd_move_f(cr->dd, &forceOutMtsLevel1->forceWithShiftForces(), wcycle);
2135 // With both nonbonded and PME offloaded a GPU on the same rank, we use
2136 // an alternating wait/reduction scheme.
2137 bool alternateGpuWait =
2138 (!c_disableAlternatingWait && stepWork.haveGpuPmeOnThisRank
2139 && simulationWork.useGpuNonbonded && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
2140 if (alternateGpuWait)
2142 alternatePmeNbGpuWaitReduce(fr->nbv.get(),
2147 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
2152 if (!alternateGpuWait && stepWork.haveGpuPmeOnThisRank)
2154 pme_gpu_wait_and_reduce(fr->pmedata,
2157 &forceOutMtsLevel1->forceWithVirial(),
2159 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]);
2162 /* Wait for local GPU NB outputs on the non-alternating wait path */
2163 if (!alternateGpuWait && stepWork.computeNonbondedForces && simulationWork.useGpuNonbonded)
2165 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
2166 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
2167 * but even with a step of 0.1 ms the difference is less than 1%
2170 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
2171 const float waitCycles = Nbnxm::gpu_wait_finish_task(
2174 AtomLocality::Local,
2175 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
2176 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
2177 forceOutNonbonded->forceWithShiftForces().shiftForces(),
2180 if (ddBalanceRegionHandler.useBalancingRegion())
2182 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
2183 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
2185 /* We measured few cycles, it could be that the kernel
2186 * and transfer finished earlier and there was no actual
2187 * wait time, only API call overhead.
2188 * Then the actual time could be anywhere between 0 and
2189 * cycles_wait_est. We will use half of cycles_wait_est.
2191 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
2193 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
2197 if (fr->nbv->emulateGpu())
2199 // NOTE: emulation kernel is not included in the balancing region,
2200 // but emulation mode does not target performance anyway
2201 wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
2206 InteractionLocality::Local,
2207 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
2211 wallcycle_stop(wcycle, WallCycleCounter::Force);
2214 // If on GPU PME-PP comms path, receive forces from PME before GPU buffer ops
2215 // TODO refactor this and unify with below default-path call to the same function
2216 if (PAR(cr) && simulationWork.haveSeparatePmeRank && simulationWork.useGpuPmePpCommunication
2217 && stepWork.computeSlowForces)
2219 /* In case of node-splitting, the PP nodes receive the long-range
2220 * forces, virial and energy from the PME nodes here.
2222 pme_receive_force_ener(fr,
2224 &forceOutMtsLevel1->forceWithVirial(),
2226 simulationWork.useGpuPmePpCommunication,
2227 stepWork.useGpuPmeFReduction,
2232 /* Do the nonbonded GPU (or emulation) force buffer reduction
2233 * on the non-alternating path. */
2234 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
2235 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
2236 if (useOrEmulateGpuNb && !alternateGpuWait)
2238 if (stepWork.useGpuFBufferOps)
2240 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2242 // TODO: move these steps as early as possible:
2243 // - CPU f H2D should be as soon as all CPU-side forces are done
2244 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
2245 // before the next CPU task that consumes the forces: vsite spread or update)
2246 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
2247 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
2248 // These should be unified.
2249 if (domainWork.haveLocalForceContribInCpuBuffer && !stepWork.useGpuFHalo)
2251 // Note: AtomLocality::All is used for the non-DD case because, as in this
2252 // case copyForcesToGpu() uses a separate stream, it allows overlap of
2253 // CPU force H2D with GPU force tasks on all streams including those in the
2254 // local stream which would otherwise be implicit dependencies for the
2255 // transfer and would not overlap.
2256 auto locality = simulationWork.havePpDomainDecomposition ? AtomLocality::Local
2257 : AtomLocality::All;
2259 stateGpu->copyForcesToGpu(forceWithShift, locality);
2262 if (stepWork.computeNonbondedForces)
2264 fr->gpuForceReduction[gmx::AtomLocality::Local]->execute();
2267 // Copy forces to host if they are needed for update or if virtual sites are enabled.
2268 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
2269 // TODO: When the output flags will be included in step workload, this copy can be combined with the
2270 // copy call done in sim_utils(...) for the output.
2271 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
2272 // they should not be copied in do_md(...) for the output.
2273 if (!simulationWork.useGpuUpdate
2274 || (simulationWork.useGpuUpdate && DOMAINDECOMP(cr) && simulationWork.useCpuPmePpCommunication)
2277 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
2278 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
2281 else if (stepWork.computeNonbondedForces)
2283 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2284 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
2288 launchGpuEndOfStepTasks(
2289 nbv, fr->listedForcesGpu.get(), fr->pmedata, enerd, *runScheduleWork, step, wcycle);
2291 if (DOMAINDECOMP(cr))
2293 dd_force_flop_stop(cr->dd, nrnb);
2296 const bool haveCombinedMtsForces = (stepWork.computeForces && simulationWork.useMts && stepWork.computeSlowForces
2297 && stepWork.combineMtsForcesBeforeHaloExchange);
2298 if (stepWork.computeForces)
2300 postProcessForceWithShiftForces(
2301 nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutMtsLevel0, vir_force, *mdatoms, *fr, vsite, stepWork);
2303 if (simulationWork.useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2305 postProcessForceWithShiftForces(
2306 nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, *mdatoms, *fr, vsite, stepWork);
2310 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
2311 if (PAR(cr) && simulationWork.haveSeparatePmeRank && simulationWork.useCpuPmePpCommunication
2312 && stepWork.computeSlowForces)
2314 /* In case of node-splitting, the PP nodes receive the long-range
2315 * forces, virial and energy from the PME nodes here.
2317 pme_receive_force_ener(fr,
2319 &forceOutMtsLevel1->forceWithVirial(),
2321 simulationWork.useGpuPmePpCommunication,
2326 if (stepWork.computeForces)
2328 /* If we don't use MTS or if we already combined the MTS forces before, we only
2329 * need to post-process one ForceOutputs object here, called forceOutCombined,
2330 * otherwise we have to post-process two outputs and then combine them.
2332 ForceOutputs& forceOutCombined = (haveCombinedMtsForces ? forceOutMts.value() : forceOutMtsLevel0);
2334 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutCombined, vir_force, mdatoms, fr, vsite, stepWork);
2336 if (simulationWork.useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2339 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, mdatoms, fr, vsite, stepWork);
2341 combineMtsForces(mdatoms->homenr,
2342 force.unpaddedArrayRef(),
2343 forceView->forceMtsCombined(),
2344 inputrec.mtsLevels[1].stepFactor);
2348 if (stepWork.computeEnergy)
2350 /* Compute the final potential energy terms */
2351 accumulatePotentialEnergies(enerd, lambda, inputrec.fepvals.get());
2353 if (!EI_TPI(inputrec.eI))
2355 checkPotentialEnergyValidity(step, *enerd, inputrec);
2359 /* In case we don't have constraints and are using GPUs, the next balancing
2360 * region starts here.
2361 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
2362 * virial calculation and COM pulling, is not thus not included in
2363 * the balance timing, which is ok as most tasks do communication.
2365 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);