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 bool havePPDomainDecomposition,
909 const SimulationWorkload& simulationWork,
910 const StepWorkload& stepWork)
912 DomainLifetimeWorkload domainWork;
913 // Note that haveSpecialForces is constant over the whole run
914 domainWork.haveSpecialForces =
915 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
916 domainWork.haveCpuListedForceWork = false;
917 domainWork.haveCpuBondedWork = false;
918 for (const auto& listedForces : fr.listedForces)
920 if (listedForces.haveCpuListedForces(*fr.fcdata))
922 domainWork.haveCpuListedForceWork = true;
924 if (listedForces.haveCpuBondeds())
926 domainWork.haveCpuBondedWork = true;
929 domainWork.haveGpuBondedWork =
930 ((fr.listedForcesGpu != nullptr) && fr.listedForcesGpu->haveInteractions());
931 // Note that haveFreeEnergyWork is constant over the whole run
932 domainWork.haveFreeEnergyWork =
933 (fr.efep != FreeEnergyPerturbationType::No && mdatoms.nPerturbed != 0);
934 // We assume we have local force work if there are CPU
935 // force tasks including PME or nonbondeds.
936 domainWork.haveCpuLocalForceWork =
937 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
938 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
939 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
940 domainWork.haveLocalForceContribInCpuBuffer =
941 domainWork.haveCpuLocalForceWork || havePPDomainDecomposition;
942 domainWork.haveNonLocalForceContribInCpuBuffer =
943 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
948 /*! \brief Set up force flag stuct from the force bitmask.
950 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
951 * \param[in] mtsLevels The multiple time-stepping levels, either empty or 2 levels
952 * \param[in] step The current MD step
953 * \param[in] simulationWork Simulation workload description.
954 * \param[in] rankHasPmeDuty If this rank computes PME.
956 * \returns New Stepworkload description.
958 static StepWorkload setupStepWorkload(const int legacyFlags,
959 ArrayRef<const gmx::MtsLevel> mtsLevels,
961 const SimulationWorkload& simulationWork,
962 const bool rankHasPmeDuty)
964 GMX_ASSERT(mtsLevels.empty() || mtsLevels.size() == 2, "Expect 0 or 2 MTS levels");
965 const bool computeSlowForces = (mtsLevels.empty() || step % mtsLevels[1].stepFactor == 0);
968 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
969 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
970 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
971 flags.computeSlowForces = computeSlowForces;
972 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
973 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
974 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
975 flags.useOnlyMtsCombinedForceBuffer = ((legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0);
976 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
977 flags.computeNonbondedForces =
978 ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded
979 && !(simulationWork.computeNonbondedAtMtsLevel1 && !computeSlowForces);
980 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
982 if (simulationWork.useGpuBufferOps)
984 GMX_ASSERT(simulationWork.useGpuNonbonded,
985 "Can only offload buffer ops if nonbonded computation is also offloaded");
987 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
988 // on virial steps the CPU reduction path is taken
989 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
990 flags.useGpuPmeFReduction = flags.computeSlowForces && flags.useGpuFBufferOps && simulationWork.useGpuPme
991 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication);
992 flags.useGpuXHalo = simulationWork.useGpuHaloExchange;
993 flags.useGpuFHalo = simulationWork.useGpuHaloExchange && flags.useGpuFBufferOps;
994 flags.haveGpuPmeOnThisRank = simulationWork.useGpuPme && rankHasPmeDuty && flags.computeSlowForces;
1000 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
1003 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
1004 gmx::ListedForcesGpu* listedForcesGpu,
1006 gmx_enerdata_t* enerd,
1007 const gmx::MdrunScheduleWorkload& runScheduleWork,
1009 gmx_wallcycle* wcycle)
1011 if (runScheduleWork.simulationWork.useGpuNonbonded && runScheduleWork.stepWork.computeNonbondedForces)
1013 /* Launch pruning before buffer clearing because the API overhead of the
1014 * clear kernel launches can leave the GPU idle while it could be running
1017 if (nbv->isDynamicPruningStepGpu(step))
1019 nbv->dispatchPruneKernelGpu(step);
1022 /* now clear the GPU outputs while we finish the step on the CPU */
1023 wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1024 wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1025 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
1026 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1027 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1030 if (runScheduleWork.stepWork.haveGpuPmeOnThisRank)
1032 pme_gpu_reinit_computation(pmedata, wcycle);
1035 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
1037 // in principle this should be included in the DD balancing region,
1038 // but generally it is infrequent so we'll omit it for the sake of
1040 listedForcesGpu->waitAccumulateEnergyTerms(enerd);
1042 listedForcesGpu->clearEnergies();
1046 //! \brief Data structure to hold dipole-related data and staging arrays
1049 //! Dipole staging for fast summing over MPI
1050 gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
1051 //! Dipole staging for states A and B (index 0 and 1 resp.)
1052 gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
1056 static void reduceAndUpdateMuTot(DipoleData* dipoleData,
1057 const t_commrec* cr,
1058 const bool haveFreeEnergy,
1059 gmx::ArrayRef<const real> lambda,
1061 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1065 gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
1066 ddBalanceRegionHandler.reopenRegionCpu();
1068 for (int i = 0; i < 2; i++)
1070 for (int j = 0; j < DIM; j++)
1072 dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
1076 if (!haveFreeEnergy)
1078 copy_rvec(dipoleData->muStateAB[0], muTotal);
1082 for (int j = 0; j < DIM; j++)
1084 muTotal[j] = (1.0 - lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)])
1085 * dipoleData->muStateAB[0][j]
1086 + lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]
1087 * dipoleData->muStateAB[1][j];
1092 /*! \brief Combines MTS level0 and level1 force buffes into a full and MTS-combined force buffer.
1094 * \param[in] numAtoms The number of atoms to combine forces for
1095 * \param[in,out] forceMtsLevel0 Input: F_level0, output: F_level0 + F_level1
1096 * \param[in,out] forceMts Input: F_level1, output: F_level0 + mtsFactor * F_level1
1097 * \param[in] mtsFactor The factor between the level0 and level1 time step
1099 static void combineMtsForces(const int numAtoms,
1100 ArrayRef<RVec> forceMtsLevel0,
1101 ArrayRef<RVec> forceMts,
1102 const real mtsFactor)
1104 const int gmx_unused numThreads = gmx_omp_nthreads_get(ModuleMultiThread::Default);
1105 #pragma omp parallel for num_threads(numThreads) schedule(static)
1106 for (int i = 0; i < numAtoms; i++)
1108 const RVec forceMtsLevel0Tmp = forceMtsLevel0[i];
1109 forceMtsLevel0[i] += forceMts[i];
1110 forceMts[i] = forceMtsLevel0Tmp + mtsFactor * forceMts[i];
1114 /*! \brief Setup for the local and non-local GPU force reductions:
1115 * reinitialization plus the registration of forces and dependencies.
1117 * \param [in] runScheduleWork Schedule workload flag structure
1118 * \param [in] cr Communication record object
1119 * \param [in] fr Force record object
1121 static void setupGpuForceReductions(gmx::MdrunScheduleWorkload* runScheduleWork,
1122 const t_commrec* cr,
1126 nonbonded_verlet_t* nbv = fr->nbv.get();
1127 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1129 // (re-)initialize local GPU force reduction
1130 const bool accumulate =
1131 runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr);
1132 const int atomStart = 0;
1133 fr->gpuForceReduction[gmx::AtomLocality::Local]->reinit(stateGpu->getForces(),
1134 nbv->getNumAtoms(AtomLocality::Local),
1135 nbv->getGridIndices(),
1138 stateGpu->fReducedOnDevice());
1140 // register forces and add dependencies
1141 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerNbnxmForce(Nbnxm::gpu_get_f(nbv->gpu_nbv));
1143 if (runScheduleWork->simulationWork.useGpuPme
1144 && (thisRankHasDuty(cr, DUTY_PME) || runScheduleWork->simulationWork.useGpuPmePpCommunication))
1146 DeviceBuffer<gmx::RVec> forcePtr =
1147 thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1148 : // PME force buffer on same GPU
1149 fr->pmePpCommGpu->getGpuForceStagingPtr(); // buffer received from other GPU
1150 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerRvecForce(forcePtr);
1152 GpuEventSynchronizer* const pmeSynchronizer =
1153 (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1154 : // PME force buffer on same GPU
1155 fr->pmePpCommGpu->getForcesReadySynchronizer()); // buffer received from other GPU
1159 GMX_ASSERT(pmeSynchronizer != nullptr, "PME force ready cuda event should not be NULL");
1160 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(pmeSynchronizer);
1164 if (runScheduleWork->domainWork.haveCpuLocalForceWork && !runScheduleWork->simulationWork.useGpuHaloExchange)
1166 // in the DD case we use the same stream for H2D and reduction, hence no explicit dependency needed
1167 if (!havePPDomainDecomposition(cr))
1169 const bool useGpuForceBufferOps = true;
1170 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1171 stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::All, useGpuForceBufferOps));
1175 if (runScheduleWork->simulationWork.useGpuHaloExchange)
1177 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1178 cr->dd->gpuHaloExchange[0][0]->getForcesReadyOnDeviceEvent());
1181 if (havePPDomainDecomposition(cr))
1183 // (re-)initialize non-local GPU force reduction
1184 const bool accumulate = runScheduleWork->domainWork.haveCpuBondedWork
1185 || runScheduleWork->domainWork.haveFreeEnergyWork;
1186 const int atomStart = dd_numHomeAtoms(*cr->dd);
1187 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->reinit(stateGpu->getForces(),
1188 nbv->getNumAtoms(AtomLocality::NonLocal),
1189 nbv->getGridIndices(),
1193 // register forces and add dependencies
1194 // in the DD case we use the same stream for H2D and reduction, hence no explicit dependency needed
1195 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->registerNbnxmForce(
1196 Nbnxm::gpu_get_f(nbv->gpu_nbv));
1201 /*! \brief Return the number of local atoms.
1203 static int getLocalAtomCount(const gmx_domdec_t* dd, const t_mdatoms& mdatoms, bool havePPDomainDecomposition)
1205 GMX_ASSERT(!(havePPDomainDecomposition && (dd == nullptr)), "Can't have PP decomposition with dd uninitialized!");
1206 return havePPDomainDecomposition ? dd_numAtomsZones(*dd) : mdatoms.homenr;
1210 void do_force(FILE* fplog,
1211 const t_commrec* cr,
1212 const gmx_multisim_t* ms,
1213 const t_inputrec& inputrec,
1215 gmx_enfrot* enforcedRotation,
1216 gmx::ImdSession* imdSession,
1220 gmx_wallcycle* wcycle,
1221 const gmx_localtop_t* top,
1223 gmx::ArrayRefWithPadding<gmx::RVec> x,
1224 const history_t* hist,
1225 gmx::ForceBuffersView* forceView,
1227 const t_mdatoms* mdatoms,
1228 gmx_enerdata_t* enerd,
1229 gmx::ArrayRef<const real> lambda,
1231 gmx::MdrunScheduleWorkload* runScheduleWork,
1232 gmx::VirtualSitesHandler* vsite,
1237 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1239 auto force = forceView->forceWithPadding();
1240 GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1241 "The size of the force buffer should be at least the number of atoms to compute "
1244 nonbonded_verlet_t* nbv = fr->nbv.get();
1245 interaction_const_t* ic = fr->ic.get();
1247 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1249 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1251 runScheduleWork->stepWork = setupStepWorkload(
1252 legacyFlags, inputrec.mtsLevels, step, simulationWork, thisRankHasDuty(cr, DUTY_PME));
1253 const StepWorkload& stepWork = runScheduleWork->stepWork;
1255 /* At a search step we need to start the first balancing region
1256 * somewhere early inside the step after communication during domain
1257 * decomposition (and not during the previous step as usual).
1259 if (stepWork.doNeighborSearch)
1261 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1264 clear_mat(vir_force);
1266 if (fr->pbcType != PbcType::No)
1268 /* Compute shift vectors every step,
1269 * because of pressure coupling or box deformation!
1271 if (stepWork.haveDynamicBox && stepWork.stateChanged)
1273 calc_shifts(box, fr->shift_vec);
1276 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1277 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1280 put_atoms_in_box_omp(fr->pbcType,
1282 x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1283 gmx_omp_nthreads_get(ModuleMultiThread::Default));
1284 inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1288 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1290 const bool pmeSendCoordinatesFromGpu =
1291 GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1292 const bool reinitGpuPmePpComms =
1293 GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1295 auto* localXReadyOnDevice = (stepWork.haveGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1296 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1297 AtomLocality::Local, simulationWork, stepWork)
1300 // Copy coordinate from the GPU if update is on the GPU and there
1301 // are forces to be computed on the CPU, or for the computation of
1302 // virial, or if host-side data will be transferred from this task
1303 // to a remote task for halo exchange or PME-PP communication. At
1304 // search steps the current coordinates are already on the host,
1305 // hence copy is not needed.
1306 const bool haveHostPmePpComms =
1307 !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1309 GMX_ASSERT(simulationWork.useGpuHaloExchange
1310 == ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange[0].empty())),
1311 "The GPU halo exchange is active, but it has not been constructed.");
1312 const bool haveHostHaloExchangeComms =
1313 havePPDomainDecomposition(cr) && !simulationWork.useGpuHaloExchange;
1315 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1316 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1317 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1318 || haveHostPmePpComms || haveHostHaloExchangeComms || simulationWork.computeMuTot))
1320 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1321 haveCopiedXFromGpu = true;
1324 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1325 // The local coordinates can be copied right away.
1326 // NOTE: Consider moving this copy to right after they are updated and constrained,
1327 // if the later is not offloaded.
1328 if (stepWork.haveGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1330 if (stepWork.doNeighborSearch)
1332 // TODO refactor this to do_md, after partitioning.
1333 stateGpu->reinit(mdatoms->homenr,
1334 getLocalAtomCount(cr->dd, *mdatoms, havePPDomainDecomposition(cr)));
1335 if (stepWork.haveGpuPmeOnThisRank)
1337 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1338 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1341 // We need to copy coordinates when:
1342 // 1. Update is not offloaded
1343 // 2. The buffers were reinitialized on search step
1344 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1346 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1347 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1351 if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces)
1353 /* Send particle coordinates to the pme nodes */
1354 if (!pmeSendCoordinatesFromGpu && !stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1356 GMX_ASSERT(haveCopiedXFromGpu,
1357 "a wait should only be triggered if copy has been scheduled");
1358 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1361 gmx_pme_send_coordinates(fr,
1364 x.unpaddedArrayRef(),
1365 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1366 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1367 (stepWork.computeVirial || stepWork.computeEnergy),
1369 simulationWork.useGpuPmePpCommunication,
1370 reinitGpuPmePpComms,
1371 pmeSendCoordinatesFromGpu,
1372 localXReadyOnDevice,
1376 if (stepWork.haveGpuPmeOnThisRank)
1378 launchPmeGpuSpread(fr->pmedata,
1381 localXReadyOnDevice,
1382 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1386 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1388 /* do gridding for pair search */
1389 if (stepWork.doNeighborSearch)
1391 if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1393 fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1396 wallcycle_start(wcycle, WallCycleCounter::NS);
1397 if (!DOMAINDECOMP(cr))
1399 const rvec vzero = { 0.0_real, 0.0_real, 0.0_real };
1400 const rvec boxDiagonal = { box[XX][XX], box[YY][YY], box[ZZ][ZZ] };
1401 wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSGridLocal);
1402 nbnxn_put_on_grid(nbv,
1408 { 0, mdatoms->homenr },
1411 x.unpaddedArrayRef(),
1414 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSGridLocal);
1418 wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSGridNonLocal);
1419 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->atomInfo, x.unpaddedArrayRef());
1420 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSGridNonLocal);
1423 nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1424 gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1427 wallcycle_stop(wcycle, WallCycleCounter::NS);
1429 /* initialize the GPU nbnxm atom data and bonded data structures */
1430 if (simulationWork.useGpuNonbonded)
1432 // Note: cycle counting only nononbondeds, GPU listed forces counts internally
1433 wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1434 wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1435 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1436 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1437 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1439 if (fr->listedForcesGpu)
1441 /* Now we put all atoms on the grid, we can assign bonded
1442 * interactions to the GPU, where the grid order is
1443 * needed. Also the xq, f and fshift device buffers have
1444 * been reallocated if needed, so the bonded code can
1445 * learn about them. */
1446 // TODO the xq, f, and fshift buffers are now shared
1447 // resources, so they should be maintained by a
1448 // higher-level object than the nb module.
1449 fr->listedForcesGpu->updateInteractionListsAndDeviceBuffers(
1450 nbv->getGridIndices(),
1452 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1453 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1454 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1458 // Need to run after the GPU-offload bonded interaction lists
1459 // are set up to be able to determine whether there is bonded work.
1460 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1461 inputrec, *fr, pull_work, ed, *mdatoms, havePPDomainDecomposition(cr), simulationWork, stepWork);
1463 wallcycle_start_nocount(wcycle, WallCycleCounter::NS);
1464 wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSSearchLocal);
1465 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1466 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1468 nbv->setupGpuShortRangeWork(fr->listedForcesGpu, InteractionLocality::Local);
1470 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSSearchLocal);
1471 wallcycle_stop(wcycle, WallCycleCounter::NS);
1473 if (stepWork.useGpuXBufferOps)
1475 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1478 if (simulationWork.useGpuBufferOps)
1480 setupGpuForceReductions(runScheduleWork, cr, fr);
1483 else if (!EI_TPI(inputrec.eI) && stepWork.computeNonbondedForces)
1485 if (stepWork.useGpuXBufferOps)
1487 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1488 nbv->convertCoordinatesGpu(AtomLocality::Local, stateGpu->getCoordinates(), localXReadyOnDevice);
1492 if (simulationWork.useGpuUpdate)
1494 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1495 GMX_ASSERT(haveCopiedXFromGpu,
1496 "a wait should only be triggered if copy has been scheduled");
1497 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1499 nbv->convertCoordinates(AtomLocality::Local, x.unpaddedArrayRef());
1503 if (simulationWork.useGpuNonbonded && (stepWork.computeNonbondedForces || domainWork.haveGpuBondedWork))
1505 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1507 wallcycle_start(wcycle, WallCycleCounter::LaunchGpu);
1508 wallcycle_sub_start(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1509 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1510 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1512 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1514 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1515 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1516 // with X buffer ops offloaded to the GPU on all but the search steps
1518 // bonded work not split into separate local and non-local, so with DD
1519 // we can only launch the kernel after non-local coordinates have been received.
1520 if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1522 fr->listedForcesGpu->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1525 /* launch local nonbonded work on GPU */
1526 wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1527 wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1528 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1529 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1530 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1533 if (stepWork.haveGpuPmeOnThisRank)
1535 // In PME GPU and mixed mode we launch FFT / gather after the
1536 // X copy/transform to allow overlap as well as after the GPU NB
1537 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1538 // the nonbonded kernel.
1539 launchPmeGpuFftAndGather(fr->pmedata,
1540 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1545 /* Communicate coordinates and sum dipole if necessary +
1546 do non-local pair search */
1547 if (havePPDomainDecomposition(cr))
1549 if (stepWork.doNeighborSearch)
1551 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1552 wallcycle_start_nocount(wcycle, WallCycleCounter::NS);
1553 wallcycle_sub_start(wcycle, WallCycleSubCounter::NBSSearchNonLocal);
1554 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1555 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1557 nbv->setupGpuShortRangeWork(fr->listedForcesGpu, InteractionLocality::NonLocal);
1558 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSSearchNonLocal);
1559 wallcycle_stop(wcycle, WallCycleCounter::NS);
1560 // TODO refactor this GPU halo exchange re-initialisation
1561 // to location in do_md where GPU halo exchange is
1562 // constructed at partitioning, after above stateGpu
1563 // re-initialization has similarly been refactored
1564 if (simulationWork.useGpuHaloExchange)
1566 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1571 if (stepWork.useGpuXHalo)
1573 // The following must be called after local setCoordinates (which records an event
1574 // when the coordinate data has been copied to the device).
1575 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1577 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1579 // non-local part of coordinate buffer must be copied back to host for CPU work
1580 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1585 if (simulationWork.useGpuUpdate)
1587 GMX_ASSERT(haveCopiedXFromGpu,
1588 "a wait should only be triggered if copy has been scheduled");
1589 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1591 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1594 if (stepWork.useGpuXBufferOps)
1596 if (!stepWork.haveGpuPmeOnThisRank && !stepWork.useGpuXHalo)
1598 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1600 nbv->convertCoordinatesGpu(AtomLocality::NonLocal,
1601 stateGpu->getCoordinates(),
1602 stateGpu->getCoordinatesReadyOnDeviceEvent(
1603 AtomLocality::NonLocal, simulationWork, stepWork));
1607 nbv->convertCoordinates(AtomLocality::NonLocal, x.unpaddedArrayRef());
1611 if (simulationWork.useGpuNonbonded)
1614 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1616 wallcycle_start(wcycle, WallCycleCounter::LaunchGpu);
1617 wallcycle_sub_start(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1618 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1619 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1620 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1623 if (domainWork.haveGpuBondedWork)
1625 fr->listedForcesGpu->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1628 /* launch non-local nonbonded tasks on GPU */
1629 wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1630 wallcycle_sub_start(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1631 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1632 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1633 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1637 if (simulationWork.useGpuNonbonded && stepWork.computeNonbondedForces)
1639 /* launch D2H copy-back F */
1640 wallcycle_start_nocount(wcycle, WallCycleCounter::LaunchGpu);
1641 wallcycle_sub_start_nocount(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1643 if (havePPDomainDecomposition(cr))
1645 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1647 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1648 wallcycle_sub_stop(wcycle, WallCycleSubCounter::LaunchGpuNonBonded);
1650 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1652 fr->listedForcesGpu->launchEnergyTransfer();
1654 wallcycle_stop(wcycle, WallCycleCounter::LaunchGpu);
1657 gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1658 if (fr->wholeMoleculeTransform)
1660 xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1663 DipoleData dipoleData;
1665 if (simulationWork.computeMuTot)
1667 const int start = 0;
1669 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch)
1671 GMX_ASSERT(haveCopiedXFromGpu,
1672 "a wait should only be triggered if copy has been scheduled");
1673 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1676 /* Calculate total (local) dipole moment in a temporary common array.
1677 * This makes it possible to sum them over nodes faster.
1679 gmx::ArrayRef<const gmx::RVec> xRef =
1680 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1684 mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
1685 : gmx::ArrayRef<real>{},
1686 mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
1687 : gmx::ArrayRef<real>{},
1688 mdatoms->nChargePerturbed != 0,
1689 dipoleData.muStaging[0],
1690 dipoleData.muStaging[1]);
1692 reduceAndUpdateMuTot(
1693 &dipoleData, cr, (fr->efep != FreeEnergyPerturbationType::No), lambda, muTotal, ddBalanceRegionHandler);
1696 /* Reset energies */
1697 reset_enerdata(enerd);
1699 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1701 wallcycle_start(wcycle, WallCycleCounter::PpDuringPme);
1702 dd_force_flop_start(cr->dd, nrnb);
1705 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1706 // this wait ensures that the D2H transfer is complete.
1707 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1708 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1710 GMX_ASSERT(haveCopiedXFromGpu, "a wait should only be triggered if copy has been scheduled");
1711 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1716 wallcycle_start(wcycle, WallCycleCounter::Rot);
1717 do_rotation(cr, enforcedRotation, box, x.unpaddedConstArrayRef(), t, step, stepWork.doNeighborSearch);
1718 wallcycle_stop(wcycle, WallCycleCounter::Rot);
1721 /* Start the force cycle counter.
1722 * Note that a different counter is used for dynamic load balancing.
1724 wallcycle_start(wcycle, WallCycleCounter::Force);
1726 /* Set up and clear force outputs:
1727 * forceOutMtsLevel0: everything except what is in the other two outputs
1728 * forceOutMtsLevel1: PME-mesh and listed-forces group 1
1729 * forceOutNonbonded: non-bonded forces
1730 * Without multiple time stepping all point to the same object.
1731 * With multiple time-stepping the use is different for MTS fast (level0 only) and slow steps.
1733 ForceOutputs forceOutMtsLevel0 = setupForceOutputs(
1734 &fr->forceHelperBuffers[0], force, domainWork, stepWork, havePPDomainDecomposition(cr), wcycle);
1736 // Force output for MTS combined forces, only set at level1 MTS steps
1737 std::optional<ForceOutputs> forceOutMts =
1738 (fr->useMts && stepWork.computeSlowForces)
1739 ? std::optional(setupForceOutputs(&fr->forceHelperBuffers[1],
1740 forceView->forceMtsCombinedWithPadding(),
1743 havePPDomainDecomposition(cr),
1747 ForceOutputs* forceOutMtsLevel1 =
1748 fr->useMts ? (stepWork.computeSlowForces ? &forceOutMts.value() : nullptr) : &forceOutMtsLevel0;
1750 const bool nonbondedAtMtsLevel1 = runScheduleWork->simulationWork.computeNonbondedAtMtsLevel1;
1752 ForceOutputs* forceOutNonbonded = nonbondedAtMtsLevel1 ? forceOutMtsLevel1 : &forceOutMtsLevel0;
1754 if (inputrec.bPull && pull_have_constraint(*pull_work))
1756 clear_pull_forces(pull_work);
1759 /* We calculate the non-bonded forces, when done on the CPU, here.
1760 * We do this before calling do_force_lowlevel, because in that
1761 * function, the listed forces are calculated before PME, which
1762 * does communication. With this order, non-bonded and listed
1763 * force calculation imbalance can be balanced out by the domain
1764 * decomposition load balancing.
1767 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1769 if (!useOrEmulateGpuNb)
1771 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1774 if (fr->efep != FreeEnergyPerturbationType::No && stepWork.computeNonbondedForces)
1776 /* Calculate the local and non-local free energy interactions here.
1777 * Happens here on the CPU both with and without GPU.
1779 nbv->dispatchFreeEnergyKernel(
1780 InteractionLocality::Local,
1781 x.unpaddedArrayRef(),
1782 &forceOutNonbonded->forceWithShiftForces(),
1783 fr->use_simd_kernels,
1790 mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
1791 : gmx::ArrayRef<real>{},
1792 mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
1793 : gmx::ArrayRef<real>{},
1794 mdatoms->typeA ? gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr)
1795 : gmx::ArrayRef<int>{},
1796 mdatoms->typeB ? gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr)
1797 : gmx::ArrayRef<int>{},
1798 inputrec.fepvals.get(),
1804 if (havePPDomainDecomposition(cr))
1806 nbv->dispatchFreeEnergyKernel(
1807 InteractionLocality::NonLocal,
1808 x.unpaddedArrayRef(),
1809 &forceOutNonbonded->forceWithShiftForces(),
1810 fr->use_simd_kernels,
1817 mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
1818 : gmx::ArrayRef<real>{},
1819 mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
1820 : gmx::ArrayRef<real>{},
1821 mdatoms->typeA ? gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr)
1822 : gmx::ArrayRef<int>{},
1823 mdatoms->typeB ? gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr)
1824 : gmx::ArrayRef<int>{},
1825 inputrec.fepvals.get(),
1833 if (stepWork.computeNonbondedForces && !useOrEmulateGpuNb)
1835 if (havePPDomainDecomposition(cr))
1837 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1840 if (stepWork.computeForces)
1842 /* Add all the non-bonded force to the normal force array.
1843 * This can be split into a local and a non-local part when overlapping
1844 * communication with calculation with domain decomposition.
1846 wallcycle_stop(wcycle, WallCycleCounter::Force);
1847 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All,
1848 forceOutNonbonded->forceWithShiftForces().force());
1849 wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
1852 /* If there are multiple fshift output buffers we need to reduce them */
1853 if (stepWork.computeVirial)
1855 /* This is not in a subcounter because it takes a
1856 negligible and constant-sized amount of time */
1857 nbnxn_atomdata_add_nbat_fshift_to_fshift(
1858 *nbv->nbat, forceOutNonbonded->forceWithShiftForces().shiftForces());
1862 // TODO Force flags should include haveFreeEnergyWork for this domain
1863 if (stepWork.useGpuXHalo && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1865 wallcycle_stop(wcycle, WallCycleCounter::Force);
1866 /* Wait for non-local coordinate data to be copied from device */
1867 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1868 wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
1871 // Compute wall interactions, when present.
1872 // Note: should be moved to special forces.
1873 if (inputrec.nwall && stepWork.computeNonbondedForces)
1875 /* foreign lambda component for walls */
1876 real dvdl_walls = do_walls(inputrec,
1879 mdatoms->typeA ? gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr)
1880 : gmx::ArrayRef<int>{},
1881 mdatoms->typeB ? gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr)
1882 : gmx::ArrayRef<int>{},
1883 mdatoms->cENER ? gmx::arrayRefFromArray(mdatoms->cENER, mdatoms->nr)
1884 : gmx::ArrayRef<unsigned short>{},
1886 mdatoms->nPerturbed,
1887 x.unpaddedConstArrayRef(),
1888 &forceOutMtsLevel0.forceWithVirial(),
1889 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1890 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR],
1892 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += dvdl_walls;
1895 if (stepWork.computeListedForces)
1897 /* Check whether we need to take into account PBC in listed interactions */
1898 bool needMolPbc = false;
1899 for (const auto& listedForces : fr->listedForces)
1901 if (listedForces.haveCpuListedForces(*fr->fcdata))
1903 needMolPbc = fr->bMolPBC;
1911 /* Since all atoms are in the rectangular or triclinic unit-cell,
1912 * only single box vector shifts (2 in x) are required.
1914 set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1917 for (int mtsIndex = 0; mtsIndex < (fr->useMts && stepWork.computeSlowForces ? 2 : 1); mtsIndex++)
1919 ListedForces& listedForces = fr->listedForces[mtsIndex];
1920 ForceOutputs& forceOut = (mtsIndex == 0 ? forceOutMtsLevel0 : *forceOutMtsLevel1);
1921 listedForces.calculate(wcycle,
1923 inputrec.fepvals.get(),
1937 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr,
1942 if (stepWork.computeSlowForces)
1944 calculateLongRangeNonbondeds(fr,
1950 x.unpaddedConstArrayRef(),
1951 &forceOutMtsLevel1->forceWithVirial(),
1955 dipoleData.muStateAB,
1957 ddBalanceRegionHandler);
1960 wallcycle_stop(wcycle, WallCycleCounter::Force);
1962 // VdW dispersion correction, only computed on master rank to avoid double counting
1963 if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1965 // Calculate long range corrections to pressure and energy
1966 const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
1967 box, lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)]);
1969 if (stepWork.computeEnergy)
1971 enerd->term[F_DISPCORR] = correction.energy;
1972 enerd->term[F_DVDL_VDW] += correction.dvdl;
1973 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += correction.dvdl;
1975 if (stepWork.computeVirial)
1977 correction.correctVirial(vir_force);
1978 enerd->term[F_PDISPCORR] = correction.pressure;
1982 computeSpecialForces(fplog,
1994 x.unpaddedArrayRef(),
1998 &forceOutMtsLevel0.forceWithVirial(),
1999 forceOutMtsLevel1 ? &forceOutMtsLevel1->forceWithVirial() : nullptr,
2002 stepWork.doNeighborSearch);
2004 if (havePPDomainDecomposition(cr) && stepWork.computeForces && stepWork.useGpuFHalo
2005 && domainWork.haveCpuLocalForceWork)
2007 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(), AtomLocality::Local);
2010 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
2011 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
2012 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFHalo),
2013 "The schedule below does not allow for nonbonded MTS with GPU halo exchange");
2014 // Will store the amount of cycles spent waiting for the GPU that
2015 // will be later used in the DLB accounting.
2016 float cycles_wait_gpu = 0;
2017 if (useOrEmulateGpuNb && stepWork.computeNonbondedForces)
2019 auto& forceWithShiftForces = forceOutNonbonded->forceWithShiftForces();
2021 /* wait for non-local forces (or calculate in emulation mode) */
2022 if (havePPDomainDecomposition(cr))
2024 if (simulationWork.useGpuNonbonded)
2026 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
2029 AtomLocality::NonLocal,
2030 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
2031 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
2032 forceWithShiftForces.shiftForces(),
2037 wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
2039 fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes, step, nrnb, wcycle);
2040 wallcycle_stop(wcycle, WallCycleCounter::Force);
2043 if (stepWork.useGpuFBufferOps)
2045 if (domainWork.haveNonLocalForceContribInCpuBuffer)
2047 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2048 AtomLocality::NonLocal);
2052 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->execute();
2054 if (!stepWork.useGpuFHalo)
2056 // copy from GPU input for dd_move_f()
2057 stateGpu->copyForcesFromGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2058 AtomLocality::NonLocal);
2063 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
2066 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
2068 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
2073 /* Combining the forces for multiple time stepping before the halo exchange, when possible,
2074 * avoids an extra halo exchange (when DD is used) and post-processing step.
2076 const bool combineMtsForcesBeforeHaloExchange =
2077 (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces && stepWork.useOnlyMtsCombinedForceBuffer
2078 && !(stepWork.computeVirial || simulationWork.useGpuNonbonded || stepWork.haveGpuPmeOnThisRank));
2079 if (combineMtsForcesBeforeHaloExchange)
2081 combineMtsForces(getLocalAtomCount(cr->dd, *mdatoms, havePPDomainDecomposition(cr)),
2082 force.unpaddedArrayRef(),
2083 forceView->forceMtsCombined(),
2084 inputrec.mtsLevels[1].stepFactor);
2087 if (havePPDomainDecomposition(cr))
2089 /* We are done with the CPU compute.
2090 * We will now communicate the non-local forces.
2091 * If we use a GPU this will overlap with GPU work, so in that case
2092 * we do not close the DD force balancing region here.
2094 ddBalanceRegionHandler.closeAfterForceComputationCpu();
2096 if (stepWork.computeForces)
2099 if (stepWork.useGpuFHalo)
2101 // If there exist CPU forces, data from halo exchange should accumulate into these
2102 bool accumulateForces = domainWork.haveCpuLocalForceWork;
2103 if (!accumulateForces)
2105 // Force halo exchange will set a subset of local atoms with remote non-local data
2106 // First clear local portion of force array, so that untouched atoms are zero
2107 stateGpu->clearForcesOnGpu(AtomLocality::Local);
2109 communicateGpuHaloForces(*cr, accumulateForces);
2113 if (stepWork.useGpuFBufferOps)
2115 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
2118 // Without MTS or with MTS at slow steps with uncombined forces we need to
2119 // communicate the fast forces
2120 if (!fr->useMts || !combineMtsForcesBeforeHaloExchange)
2122 dd_move_f(cr->dd, &forceOutMtsLevel0.forceWithShiftForces(), wcycle);
2124 // With MTS we need to communicate the slow or combined (in forceOutMtsLevel1) forces
2125 if (fr->useMts && stepWork.computeSlowForces)
2127 dd_move_f(cr->dd, &forceOutMtsLevel1->forceWithShiftForces(), wcycle);
2133 // With both nonbonded and PME offloaded a GPU on the same rank, we use
2134 // an alternating wait/reduction scheme.
2135 bool alternateGpuWait =
2136 (!c_disableAlternatingWait && stepWork.haveGpuPmeOnThisRank
2137 && simulationWork.useGpuNonbonded && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
2138 if (alternateGpuWait)
2140 alternatePmeNbGpuWaitReduce(fr->nbv.get(),
2145 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
2150 if (!alternateGpuWait && stepWork.haveGpuPmeOnThisRank)
2152 pme_gpu_wait_and_reduce(fr->pmedata,
2155 &forceOutMtsLevel1->forceWithVirial(),
2157 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]);
2160 /* Wait for local GPU NB outputs on the non-alternating wait path */
2161 if (!alternateGpuWait && stepWork.computeNonbondedForces && simulationWork.useGpuNonbonded)
2163 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
2164 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
2165 * but even with a step of 0.1 ms the difference is less than 1%
2168 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
2169 const float waitCycles = Nbnxm::gpu_wait_finish_task(
2172 AtomLocality::Local,
2173 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
2174 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
2175 forceOutNonbonded->forceWithShiftForces().shiftForces(),
2178 if (ddBalanceRegionHandler.useBalancingRegion())
2180 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
2181 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
2183 /* We measured few cycles, it could be that the kernel
2184 * and transfer finished earlier and there was no actual
2185 * wait time, only API call overhead.
2186 * Then the actual time could be anywhere between 0 and
2187 * cycles_wait_est. We will use half of cycles_wait_est.
2189 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
2191 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
2195 if (fr->nbv->emulateGpu())
2197 // NOTE: emulation kernel is not included in the balancing region,
2198 // but emulation mode does not target performance anyway
2199 wallcycle_start_nocount(wcycle, WallCycleCounter::Force);
2204 InteractionLocality::Local,
2205 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
2209 wallcycle_stop(wcycle, WallCycleCounter::Force);
2212 // If on GPU PME-PP comms path, receive forces from PME before GPU buffer ops
2213 // TODO refactor this and unify with below default-path call to the same function
2214 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces
2215 && simulationWork.useGpuPmePpCommunication)
2217 /* In case of node-splitting, the PP nodes receive the long-range
2218 * forces, virial and energy from the PME nodes here.
2220 pme_receive_force_ener(fr,
2222 &forceOutMtsLevel1->forceWithVirial(),
2224 simulationWork.useGpuPmePpCommunication,
2225 stepWork.useGpuPmeFReduction,
2230 /* Do the nonbonded GPU (or emulation) force buffer reduction
2231 * on the non-alternating path. */
2232 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
2233 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
2234 if (useOrEmulateGpuNb && !alternateGpuWait)
2236 if (stepWork.useGpuFBufferOps)
2238 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2240 // TODO: move these steps as early as possible:
2241 // - CPU f H2D should be as soon as all CPU-side forces are done
2242 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
2243 // before the next CPU task that consumes the forces: vsite spread or update)
2244 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
2245 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
2246 // These should be unified.
2247 if (domainWork.haveLocalForceContribInCpuBuffer && !stepWork.useGpuFHalo)
2249 // Note: AtomLocality::All is used for the non-DD case because, as in this
2250 // case copyForcesToGpu() uses a separate stream, it allows overlap of
2251 // CPU force H2D with GPU force tasks on all streams including those in the
2252 // local stream which would otherwise be implicit dependencies for the
2253 // transfer and would not overlap.
2254 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
2256 stateGpu->copyForcesToGpu(forceWithShift, locality);
2259 if (stepWork.computeNonbondedForces)
2261 fr->gpuForceReduction[gmx::AtomLocality::Local]->execute();
2264 // Copy forces to host if they are needed for update or if virtual sites are enabled.
2265 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
2266 // TODO: When the output flags will be included in step workload, this copy can be combined with the
2267 // copy call done in sim_utils(...) for the output.
2268 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
2269 // they should not be copied in do_md(...) for the output.
2270 if (!simulationWork.useGpuUpdate
2271 || (simulationWork.useGpuUpdate && DOMAINDECOMP(cr) && haveHostPmePpComms) || vsite)
2273 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
2274 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
2277 else if (stepWork.computeNonbondedForces)
2279 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2280 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
2284 launchGpuEndOfStepTasks(nbv, fr->listedForcesGpu, fr->pmedata, enerd, *runScheduleWork, step, wcycle);
2286 if (DOMAINDECOMP(cr))
2288 dd_force_flop_stop(cr->dd, nrnb);
2291 const bool haveCombinedMtsForces = (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2292 && combineMtsForcesBeforeHaloExchange);
2293 if (stepWork.computeForces)
2295 postProcessForceWithShiftForces(
2296 nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutMtsLevel0, vir_force, *mdatoms, *fr, vsite, stepWork);
2298 if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2300 postProcessForceWithShiftForces(
2301 nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, *mdatoms, *fr, vsite, stepWork);
2305 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
2306 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
2307 && stepWork.computeSlowForces)
2309 /* In case of node-splitting, the PP nodes receive the long-range
2310 * forces, virial and energy from the PME nodes here.
2312 pme_receive_force_ener(fr,
2314 &forceOutMtsLevel1->forceWithVirial(),
2316 simulationWork.useGpuPmePpCommunication,
2321 if (stepWork.computeForces)
2323 /* If we don't use MTS or if we already combined the MTS forces before, we only
2324 * need to post-process one ForceOutputs object here, called forceOutCombined,
2325 * otherwise we have to post-process two outputs and then combine them.
2327 ForceOutputs& forceOutCombined = (haveCombinedMtsForces ? forceOutMts.value() : forceOutMtsLevel0);
2329 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutCombined, vir_force, mdatoms, fr, vsite, stepWork);
2331 if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2334 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, mdatoms, fr, vsite, stepWork);
2336 combineMtsForces(mdatoms->homenr,
2337 force.unpaddedArrayRef(),
2338 forceView->forceMtsCombined(),
2339 inputrec.mtsLevels[1].stepFactor);
2343 if (stepWork.computeEnergy)
2345 /* Compute the final potential energy terms */
2346 accumulatePotentialEnergies(enerd, lambda, inputrec.fepvals.get());
2348 if (!EI_TPI(inputrec.eI))
2350 checkPotentialEnergyValidity(step, *enerd, inputrec);
2354 /* In case we don't have constraints and are using GPUs, the next balancing
2355 * region starts here.
2356 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
2357 * virial calculation and COM pulling, is not thus not included in
2358 * the balance timing, which is ok as most tasks do communication.
2360 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);