2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013-2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "gromacs/awh/awh.h"
49 #include "gromacs/domdec/dlbtiming.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/domdec/gpuhaloexchange.h"
53 #include "gromacs/domdec/partition.h"
54 #include "gromacs/essentialdynamics/edsam.h"
55 #include "gromacs/ewald/pme.h"
56 #include "gromacs/ewald/pme_pp.h"
57 #include "gromacs/ewald/pme_pp_comm_gpu.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
60 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
61 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
62 #include "gromacs/gpu_utils/gpu_utils.h"
63 #include "gromacs/imd/imd.h"
64 #include "gromacs/listed_forces/disre.h"
65 #include "gromacs/listed_forces/gpubonded.h"
66 #include "gromacs/listed_forces/listed_forces.h"
67 #include "gromacs/listed_forces/orires.h"
68 #include "gromacs/math/arrayrefwithpadding.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/units.h"
71 #include "gromacs/math/vec.h"
72 #include "gromacs/math/vecdump.h"
73 #include "gromacs/mdlib/calcmu.h"
74 #include "gromacs/mdlib/calcvir.h"
75 #include "gromacs/mdlib/constr.h"
76 #include "gromacs/mdlib/dispersioncorrection.h"
77 #include "gromacs/mdlib/enerdata_utils.h"
78 #include "gromacs/mdlib/force.h"
79 #include "gromacs/mdlib/force_flags.h"
80 #include "gromacs/mdlib/forcerec.h"
81 #include "gromacs/mdlib/gmx_omp_nthreads.h"
82 #include "gromacs/mdlib/update.h"
83 #include "gromacs/mdlib/vsite.h"
84 #include "gromacs/mdlib/wholemoleculetransform.h"
85 #include "gromacs/mdtypes/commrec.h"
86 #include "gromacs/mdtypes/enerdata.h"
87 #include "gromacs/mdtypes/forcebuffers.h"
88 #include "gromacs/mdtypes/forceoutput.h"
89 #include "gromacs/mdtypes/forcerec.h"
90 #include "gromacs/mdtypes/iforceprovider.h"
91 #include "gromacs/mdtypes/inputrec.h"
92 #include "gromacs/mdtypes/md_enums.h"
93 #include "gromacs/mdtypes/mdatom.h"
94 #include "gromacs/mdtypes/simulation_workload.h"
95 #include "gromacs/mdtypes/state.h"
96 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
97 #include "gromacs/nbnxm/gpu_data_mgmt.h"
98 #include "gromacs/nbnxm/nbnxm.h"
99 #include "gromacs/nbnxm/nbnxm_gpu.h"
100 #include "gromacs/pbcutil/ishift.h"
101 #include "gromacs/pbcutil/pbc.h"
102 #include "gromacs/pulling/pull.h"
103 #include "gromacs/pulling/pull_rotation.h"
104 #include "gromacs/timing/cyclecounter.h"
105 #include "gromacs/timing/gpu_timing.h"
106 #include "gromacs/timing/wallcycle.h"
107 #include "gromacs/timing/wallcyclereporting.h"
108 #include "gromacs/timing/walltime_accounting.h"
109 #include "gromacs/topology/topology.h"
110 #include "gromacs/utility/arrayref.h"
111 #include "gromacs/utility/basedefinitions.h"
112 #include "gromacs/utility/cstringutil.h"
113 #include "gromacs/utility/exceptions.h"
114 #include "gromacs/utility/fatalerror.h"
115 #include "gromacs/utility/fixedcapacityvector.h"
116 #include "gromacs/utility/gmxassert.h"
117 #include "gromacs/utility/gmxmpi.h"
118 #include "gromacs/utility/logger.h"
119 #include "gromacs/utility/smalloc.h"
120 #include "gromacs/utility/strconvert.h"
121 #include "gromacs/utility/sysinfo.h"
124 using gmx::AtomLocality;
125 using gmx::DomainLifetimeWorkload;
126 using gmx::ForceOutputs;
127 using gmx::ForceWithShiftForces;
128 using gmx::InteractionLocality;
130 using gmx::SimulationWorkload;
131 using gmx::StepWorkload;
133 // TODO: this environment variable allows us to verify before release
134 // that on less common architectures the total cost of polling is not larger than
135 // a blocking wait (so polling does not introduce overhead when the static
136 // PME-first ordering would suffice).
137 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
139 static void sum_forces(ArrayRef<RVec> f, ArrayRef<const RVec> forceToAdd)
141 GMX_ASSERT(f.size() >= forceToAdd.size(), "Accumulation buffer should be sufficiently large");
142 const int end = forceToAdd.size();
144 int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
145 #pragma omp parallel for num_threads(nt) schedule(static)
146 for (int i = 0; i < end; i++)
148 rvec_inc(f[i], forceToAdd[i]);
152 static void calc_virial(int start,
155 const gmx::ForceWithShiftForces& forceWithShiftForces,
159 const t_forcerec* fr,
162 /* The short-range virial from surrounding boxes */
163 const rvec* fshift = as_rvec_array(forceWithShiftForces.shiftForces().data());
164 calc_vir(SHIFTS, fr->shift_vec, fshift, vir_part, pbcType == PbcType::Screw, box);
165 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
167 /* Calculate partial virial, for local atoms only, based on short range.
168 * Total virial is computed in global_stat, called from do_md
170 const rvec* f = as_rvec_array(forceWithShiftForces.force().data());
171 f_calc_vir(start, start + homenr, x, f, vir_part, box);
172 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
176 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
180 static void pull_potential_wrapper(const t_commrec* cr,
181 const t_inputrec* ir,
183 gmx::ArrayRef<const gmx::RVec> x,
184 gmx::ForceWithVirial* force,
185 const t_mdatoms* mdatoms,
186 gmx_enerdata_t* enerd,
190 gmx_wallcycle_t wcycle)
195 /* Calculate the center of mass forces, this requires communication,
196 * which is why pull_potential is called close to other communication.
198 wallcycle_start(wcycle, ewcPULLPOT);
199 set_pbc(&pbc, ir->pbcType, box);
201 enerd->term[F_COM_PULL] +=
202 pull_potential(pull_work, mdatoms->massT, &pbc, cr, t, lambda[efptRESTRAINT],
203 as_rvec_array(x.data()), force, &dvdl);
204 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
205 wallcycle_stop(wcycle, ewcPULLPOT);
208 static void pme_receive_force_ener(t_forcerec* fr,
210 gmx::ForceWithVirial* forceWithVirial,
211 gmx_enerdata_t* enerd,
212 bool useGpuPmePpComms,
213 bool receivePmeForceToGpu,
214 gmx_wallcycle_t wcycle)
216 real e_q, e_lj, dvdl_q, dvdl_lj;
217 float cycles_ppdpme, cycles_seppme;
219 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
220 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
222 /* In case of node-splitting, the PP nodes receive the long-range
223 * forces, virial and energy from the PME nodes here.
225 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
228 gmx_pme_receive_f(fr->pmePpCommGpu.get(), cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
229 useGpuPmePpComms, receivePmeForceToGpu, &cycles_seppme);
230 enerd->term[F_COUL_RECIP] += e_q;
231 enerd->term[F_LJ_RECIP] += e_lj;
232 enerd->dvdl_lin[efptCOUL] += dvdl_q;
233 enerd->dvdl_lin[efptVDW] += dvdl_lj;
237 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
239 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
242 static void print_large_forces(FILE* fp,
247 ArrayRef<const RVec> x,
248 ArrayRef<const RVec> f)
250 real force2Tolerance = gmx::square(forceTolerance);
251 gmx::index numNonFinite = 0;
252 for (int i = 0; i < md->homenr; i++)
254 real force2 = norm2(f[i]);
255 bool nonFinite = !std::isfinite(force2);
256 if (force2 >= force2Tolerance || nonFinite)
258 fprintf(fp, "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n", step,
259 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
266 if (numNonFinite > 0)
268 /* Note that with MPI this fatal call on one rank might interrupt
269 * the printing on other ranks. But we can only avoid that with
270 * an expensive MPI barrier that we would need at each step.
272 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
276 //! When necessary, spreads forces on vsites and computes the virial for \p forceOutputs->forceWithShiftForces()
277 static void postProcessForceWithShiftForces(t_nrnb* nrnb,
278 gmx_wallcycle_t wcycle,
280 ArrayRef<const RVec> x,
281 ForceOutputs* forceOutputs,
283 const t_mdatoms& mdatoms,
284 const t_forcerec& fr,
285 gmx::VirtualSitesHandler* vsite,
286 const StepWorkload& stepWork)
288 ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
290 /* If we have NoVirSum forces, but we do not calculate the virial,
291 * we later sum the forceWithShiftForces buffer together with
292 * the noVirSum buffer and spread the combined vsite forces at once.
294 if (vsite && (!forceOutputs->haveForceWithVirial() || stepWork.computeVirial))
296 using VirialHandling = gmx::VirtualSitesHandler::VirialHandling;
298 auto f = forceWithShiftForces.force();
299 auto fshift = forceWithShiftForces.shiftForces();
300 const VirialHandling virialHandling =
301 (stepWork.computeVirial ? VirialHandling::Pbc : VirialHandling::None);
302 vsite->spreadForces(x, f, virialHandling, fshift, nullptr, nrnb, box, wcycle);
303 forceWithShiftForces.haveSpreadVsiteForces() = true;
306 if (stepWork.computeVirial)
308 /* Calculation of the virial must be done after vsites! */
309 calc_virial(0, mdatoms.homenr, as_rvec_array(x.data()), forceWithShiftForces, vir_force,
310 box, nrnb, &fr, fr.pbcType);
314 //! Spread, compute virial for and sum forces, when necessary
315 static void postProcessForces(const t_commrec* cr,
318 gmx_wallcycle_t wcycle,
320 ArrayRef<const RVec> x,
321 ForceOutputs* forceOutputs,
323 const t_mdatoms* mdatoms,
324 const t_forcerec* fr,
325 gmx::VirtualSitesHandler* vsite,
326 const StepWorkload& stepWork)
328 // Extract the final output force buffer, which is also the buffer for forces with shift forces
329 ArrayRef<RVec> f = forceOutputs->forceWithShiftForces().force();
331 if (forceOutputs->haveForceWithVirial())
333 auto& forceWithVirial = forceOutputs->forceWithVirial();
337 /* Spread the mesh force on virtual sites to the other particles...
338 * This is parallellized. MPI communication is performed
339 * if the constructing atoms aren't local.
341 GMX_ASSERT(!stepWork.computeVirial || f.data() != forceWithVirial.force_.data(),
342 "We need separate force buffers for shift and virial forces when "
343 "computing the virial");
344 GMX_ASSERT(!stepWork.computeVirial
345 || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
346 "We should spread the force with shift forces separately when computing "
348 const gmx::VirtualSitesHandler::VirialHandling virialHandling =
349 (stepWork.computeVirial ? gmx::VirtualSitesHandler::VirialHandling::NonLinear
350 : gmx::VirtualSitesHandler::VirialHandling::None);
351 matrix virial = { { 0 } };
352 vsite->spreadForces(x, forceWithVirial.force_, virialHandling, {}, virial, nrnb, box, wcycle);
353 forceWithVirial.addVirialContribution(virial);
356 if (stepWork.computeVirial)
358 /* Now add the forces, this is local */
359 sum_forces(f, forceWithVirial.force_);
361 /* Add the direct virial contributions */
363 forceWithVirial.computeVirial_,
364 "forceWithVirial should request virial computation when we request the virial");
365 m_add(vir_force, forceWithVirial.getVirial(), vir_force);
369 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
375 GMX_ASSERT(vsite == nullptr || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
376 "We should have spread the vsite forces (earlier)");
379 if (fr->print_force >= 0)
381 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
385 static void do_nb_verlet(t_forcerec* fr,
386 const interaction_const_t* ic,
387 gmx_enerdata_t* enerd,
388 const StepWorkload& stepWork,
389 const InteractionLocality ilocality,
393 gmx_wallcycle_t wcycle)
395 if (!stepWork.computeNonbondedForces)
397 /* skip non-bonded calculation */
401 nonbonded_verlet_t* nbv = fr->nbv.get();
403 /* GPU kernel launch overhead is already timed separately */
404 if (fr->cutoff_scheme != ecutsVERLET)
406 gmx_incons("Invalid cut-off scheme passed!");
411 /* When dynamic pair-list pruning is requested, we need to prune
412 * at nstlistPrune steps.
414 if (nbv->isDynamicPruningStepCpu(step))
416 /* Prune the pair-list beyond fr->ic->rlistPrune using
417 * the current coordinates of the atoms.
419 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
420 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
421 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
425 nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
428 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
430 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, v.ssize());
432 /* Note that we would like to avoid this conditional by putting it
433 * into the omp pragma instead, but then we still take the full
434 * omp parallel for overhead (at least with gcc5).
436 if (!useOpenmpThreading || nth == 1)
445 #pragma omp parallel for num_threads(nth) schedule(static)
446 for (gmx::index i = 0; i < v.ssize(); i++)
453 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
455 * \param groupOptions Group options, containing T-coupling options
457 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
459 real nrdfCoupled = 0;
460 real nrdfUncoupled = 0;
461 real kineticEnergy = 0;
462 for (int g = 0; g < groupOptions.ngtc; g++)
464 if (groupOptions.tau_t[g] >= 0)
466 nrdfCoupled += groupOptions.nrdf[g];
467 kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * BOLTZ;
471 nrdfUncoupled += groupOptions.nrdf[g];
475 /* This conditional with > also catches nrdf=0 */
476 if (nrdfCoupled > nrdfUncoupled)
478 return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
486 /*! \brief This routine checks that the potential energy is finite.
488 * Always checks that the potential energy is finite. If step equals
489 * inputrec.init_step also checks that the magnitude of the potential energy
490 * is reasonable. Terminates with a fatal error when a check fails.
491 * Note that passing this check does not guarantee finite forces,
492 * since those use slightly different arithmetics. But in most cases
493 * there is just a narrow coordinate range where forces are not finite
494 * and energies are finite.
496 * \param[in] step The step number, used for checking and printing
497 * \param[in] enerd The energy data; the non-bonded group energies need to be added to
498 * enerd.term[F_EPOT] before calling this routine \param[in] inputrec The input record
500 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
502 /* Threshold valid for comparing absolute potential energy against
503 * the kinetic energy. Normally one should not consider absolute
504 * potential energy values, but with a factor of one million
505 * we should never get false positives.
507 constexpr real c_thresholdFactor = 1e6;
509 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
510 real averageKineticEnergy = 0;
511 /* We only check for large potential energy at the initial step,
512 * because that is by far the most likely step for this too occur
513 * and because computing the average kinetic energy is not free.
514 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
515 * before they become NaN.
517 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
519 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
522 if (energyIsNotFinite
523 || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
528 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
529 "contributions to the energy are %g and %g, respectively. A %s potential energy "
530 "can be caused by overlapping interactions in bonded interactions or very large%s "
531 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
532 "configuration, incorrect interactions or parameters in the topology.",
533 step, enerd.term[F_EPOT], energyIsNotFinite ? "not finite" : "extremely high",
534 enerd.term[F_LJ], enerd.term[F_COUL_SR],
535 energyIsNotFinite ? "non-finite" : "very high", energyIsNotFinite ? " or Nan" : "");
539 /*! \brief Return true if there are special forces computed this step.
541 * The conditionals exactly correspond to those in computeSpecialForces().
543 static bool haveSpecialForces(const t_inputrec& inputrec,
544 const gmx::ForceProviders& forceProviders,
545 const pull_t* pull_work,
546 const bool computeForces,
550 return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
551 (inputrec.bPull && pull_have_potential(pull_work)) || // pull
552 inputrec.bRot || // enforced rotation
553 (ed != nullptr) || // flooding
554 (inputrec.bIMD && computeForces)); // IMD
557 /*! \brief Compute forces and/or energies for special algorithms
559 * The intention is to collect all calls to algorithms that compute
560 * forces on local atoms only and that do not contribute to the local
561 * virial sum (but add their virial contribution separately).
562 * Eventually these should likely all become ForceProviders.
563 * Within this function the intention is to have algorithms that do
564 * global communication at the end, so global barriers within the MD loop
565 * are as close together as possible.
567 * \param[in] fplog The log file
568 * \param[in] cr The communication record
569 * \param[in] inputrec The input record
570 * \param[in] awh The Awh module (nullptr if none in use).
571 * \param[in] enforcedRotation Enforced rotation module.
572 * \param[in] imdSession The IMD session
573 * \param[in] pull_work The pull work structure.
574 * \param[in] step The current MD step
575 * \param[in] t The current time
576 * \param[in,out] wcycle Wallcycle accounting struct
577 * \param[in,out] forceProviders Pointer to a list of force providers
578 * \param[in] box The unit cell
579 * \param[in] x The coordinates
580 * \param[in] mdatoms Per atom properties
581 * \param[in] lambda Array of free-energy lambda values
582 * \param[in] stepWork Step schedule flags
583 * \param[in,out] forceWithVirial Force and virial buffers
584 * \param[in,out] enerd Energy buffer
585 * \param[in,out] ed Essential dynamics pointer
586 * \param[in] didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
588 * \todo Remove didNeighborSearch, which is used incorrectly.
589 * \todo Convert all other algorithms called here to ForceProviders.
591 static void computeSpecialForces(FILE* fplog,
593 const t_inputrec* inputrec,
595 gmx_enfrot* enforcedRotation,
596 gmx::ImdSession* imdSession,
600 gmx_wallcycle_t wcycle,
601 gmx::ForceProviders* forceProviders,
603 gmx::ArrayRef<const gmx::RVec> x,
604 const t_mdatoms* mdatoms,
605 gmx::ArrayRef<const real> lambda,
606 const StepWorkload& stepWork,
607 gmx::ForceWithVirial* forceWithVirial,
608 gmx_enerdata_t* enerd,
610 bool didNeighborSearch)
612 /* NOTE: Currently all ForceProviders only provide forces.
613 * When they also provide energies, remove this conditional.
615 if (stepWork.computeForces)
617 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
618 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
620 /* Collect forces from modules */
621 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
624 if (inputrec->bPull && pull_have_potential(pull_work))
626 pull_potential_wrapper(cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work,
627 lambda.data(), t, wcycle);
631 const bool needForeignEnergyDifferences = awh->needForeignEnergyDifferences(step);
632 std::vector<double> foreignLambdaDeltaH, foreignLambdaDhDl;
633 if (needForeignEnergyDifferences)
635 enerd->foreignLambdaTerms.finalizePotentialContributions(enerd->dvdl_lin, lambda,
637 std::tie(foreignLambdaDeltaH, foreignLambdaDhDl) = enerd->foreignLambdaTerms.getTerms(cr);
640 enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(
641 inputrec->pbcType, mdatoms->massT, foreignLambdaDeltaH, foreignLambdaDhDl, box,
642 forceWithVirial, t, step, wcycle, fplog);
645 rvec* f = as_rvec_array(forceWithVirial->force_.data());
647 /* Add the forces from enforced rotation potentials (if any) */
650 wallcycle_start(wcycle, ewcROTadd);
651 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
652 wallcycle_stop(wcycle, ewcROTadd);
657 /* Note that since init_edsam() is called after the initialization
658 * of forcerec, edsam doesn't request the noVirSum force buffer.
659 * Thus if no other algorithm (e.g. PME) requires it, the forces
660 * here will contribute to the virial.
662 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
665 /* Add forces from interactive molecular dynamics (IMD), if any */
666 if (inputrec->bIMD && stepWork.computeForces)
668 imdSession->applyForces(f);
672 /*! \brief Launch the prepare_step and spread stages of PME GPU.
674 * \param[in] pmedata The PME structure
675 * \param[in] box The box matrix
676 * \param[in] stepWork Step schedule flags
677 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in the device memory.
678 * \param[in] lambdaQ The Coulomb lambda of the current state.
679 * \param[in] wcycle The wallcycle structure
681 static inline void launchPmeGpuSpread(gmx_pme_t* pmedata,
683 const StepWorkload& stepWork,
684 GpuEventSynchronizer* xReadyOnDevice,
686 gmx_wallcycle_t wcycle)
688 pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
689 pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle, lambdaQ);
692 /*! \brief Launch the FFT and gather stages of PME GPU
694 * This function only implements setting the output forces (no accumulation).
696 * \param[in] pmedata The PME structure
697 * \param[in] lambdaQ The Coulomb lambda of the current system state.
698 * \param[in] wcycle The wallcycle structure
699 * \param[in] stepWork Step schedule flags
701 static void launchPmeGpuFftAndGather(gmx_pme_t* pmedata,
703 gmx_wallcycle_t wcycle,
704 const gmx::StepWorkload& stepWork)
706 pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
707 pme_gpu_launch_gather(pmedata, wcycle, lambdaQ);
711 * Polling wait for either of the PME or nonbonded GPU tasks.
713 * Instead of a static order in waiting for GPU tasks, this function
714 * polls checking which of the two tasks completes first, and does the
715 * associated force buffer reduction overlapped with the other task.
716 * By doing that, unlike static scheduling order, it can always overlap
717 * one of the reductions, regardless of the GPU task completion order.
719 * \param[in] nbv Nonbonded verlet structure
720 * \param[in,out] pmedata PME module data
721 * \param[in,out] forceOutputs Output buffer for the forces and virial
722 * \param[in,out] enerd Energy data structure results are reduced into
723 * \param[in] lambdaQ The Coulomb lambda of the current system state.
724 * \param[in] stepWork Step schedule flags
725 * \param[in] wcycle The wallcycle structure
727 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
729 gmx::ForceOutputs* forceOutputs,
730 gmx_enerdata_t* enerd,
732 const StepWorkload& stepWork,
733 gmx_wallcycle_t wcycle)
735 bool isPmeGpuDone = false;
736 bool isNbGpuDone = false;
739 gmx::ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
740 gmx::ForceWithVirial& forceWithVirial = forceOutputs->forceWithVirial();
742 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
744 while (!isPmeGpuDone || !isNbGpuDone)
748 GpuTaskCompletion completionType =
749 (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
750 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, stepWork, wcycle, &forceWithVirial,
751 enerd, lambdaQ, completionType);
756 GpuTaskCompletion completionType =
757 (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
758 isNbGpuDone = Nbnxm::gpu_try_finish_task(
759 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
760 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(),
761 completionType, wcycle);
765 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShiftForces.force());
771 /*! \brief Set up the different force buffers; also does clearing.
773 * \param[in] forceHelperBuffers Helper force buffers
774 * \param[in] pull_work The pull work object.
775 * \param[in] inputrec input record
776 * \param[in] force force array
777 * \param[in] stepWork Step schedule flags
778 * \param[out] wcycle wallcycle recording structure
780 * \returns Cleared force output structure
782 static ForceOutputs setupForceOutputs(ForceHelperBuffers* forceHelperBuffers,
784 const t_inputrec& inputrec,
785 gmx::ArrayRefWithPadding<gmx::RVec> force,
786 const StepWorkload& stepWork,
787 gmx_wallcycle_t wcycle)
789 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
791 /* NOTE: We assume fr->shiftForces is all zeros here */
792 gmx::ForceWithShiftForces forceWithShiftForces(force, stepWork.computeVirial,
793 forceHelperBuffers->shiftForces());
795 if (stepWork.computeForces)
797 /* Clear the short- and long-range forces */
798 clearRVecs(forceWithShiftForces.force(), true);
800 /* Clear the shift forces */
801 clearRVecs(forceWithShiftForces.shiftForces(), false);
804 /* If we need to compute the virial, we might need a separate
805 * force buffer for algorithms for which the virial is calculated
806 * directly, such as PME. Otherwise, forceWithVirial uses the
807 * the same force (f in legacy calls) buffer as other algorithms.
809 const bool useSeparateForceWithVirialBuffer =
810 (stepWork.computeForces
811 && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
812 /* forceWithVirial uses the local atom range only */
813 gmx::ForceWithVirial forceWithVirial(
814 useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
815 : force.unpaddedArrayRef(),
816 stepWork.computeVirial);
818 if (useSeparateForceWithVirialBuffer)
820 /* TODO: update comment
821 * We only compute forces on local atoms. Note that vsites can
822 * spread to non-local atoms, but that part of the buffer is
823 * cleared separately in the vsite spreading code.
825 clearRVecs(forceWithVirial.force_, true);
828 if (inputrec.bPull && pull_have_constraint(pull_work))
830 clear_pull_forces(pull_work);
833 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
835 return ForceOutputs(forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(),
840 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
842 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec& inputrec,
843 const t_forcerec& fr,
844 const pull_t* pull_work,
846 const t_mdatoms& mdatoms,
847 const SimulationWorkload& simulationWork,
848 const StepWorkload& stepWork)
850 DomainLifetimeWorkload domainWork;
851 // Note that haveSpecialForces is constant over the whole run
852 domainWork.haveSpecialForces =
853 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
854 domainWork.haveCpuListedForceWork = false;
855 domainWork.haveCpuBondedWork = false;
856 for (const auto& listedForces : fr.listedForces)
858 if (listedForces.haveCpuListedForces(*fr.fcdata))
860 domainWork.haveCpuListedForceWork = true;
862 if (listedForces.haveCpuBondeds())
864 domainWork.haveCpuBondedWork = true;
867 domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
868 // Note that haveFreeEnergyWork is constant over the whole run
869 domainWork.haveFreeEnergyWork = (fr.efep != efepNO && mdatoms.nPerturbed != 0);
870 // We assume we have local force work if there are CPU
871 // force tasks including PME or nonbondeds.
872 domainWork.haveCpuLocalForceWork =
873 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
874 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
875 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
880 /*! \brief Set up force flag stuct from the force bitmask.
882 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
883 * \param[in] isNonbondedOn Global override, if false forces to turn off all nonbonded calculation.
884 * \param[in] simulationWork Simulation workload description.
885 * \param[in] rankHasPmeDuty If this rank computes PME.
887 * \returns New Stepworkload description.
889 static StepWorkload setupStepWorkload(const int legacyFlags,
890 const bool isNonbondedOn,
891 const SimulationWorkload& simulationWork,
892 const bool rankHasPmeDuty)
895 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
896 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
897 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
898 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
899 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
900 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
901 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
902 flags.computeNonbondedForces = ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && isNonbondedOn;
903 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
905 if (simulationWork.useGpuBufferOps)
907 GMX_ASSERT(simulationWork.useGpuNonbonded,
908 "Can only offload buffer ops if nonbonded computation is also offloaded");
910 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
911 // on virial steps the CPU reduction path is taken
912 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
913 flags.useGpuPmeFReduction = flags.useGpuFBufferOps
914 && (simulationWork.useGpuPme
915 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication));
921 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
923 * TODO: eliminate \p useGpuPmeOnThisRank when this is
924 * incorporated in DomainLifetimeWorkload.
926 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
927 gmx::GpuBonded* gpuBonded,
929 gmx_enerdata_t* enerd,
930 const gmx::MdrunScheduleWorkload& runScheduleWork,
931 bool useGpuPmeOnThisRank,
933 gmx_wallcycle_t wcycle)
935 if (runScheduleWork.simulationWork.useGpuNonbonded)
937 /* Launch pruning before buffer clearing because the API overhead of the
938 * clear kernel launches can leave the GPU idle while it could be running
941 if (nbv->isDynamicPruningStepGpu(step))
943 nbv->dispatchPruneKernelGpu(step);
946 /* now clear the GPU outputs while we finish the step on the CPU */
947 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
948 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
949 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
950 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
951 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
954 if (useGpuPmeOnThisRank)
956 pme_gpu_reinit_computation(pmedata, wcycle);
959 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
961 // in principle this should be included in the DD balancing region,
962 // but generally it is infrequent so we'll omit it for the sake of
964 gpuBonded->waitAccumulateEnergyTerms(enerd);
966 gpuBonded->clearEnergies();
970 //! \brief Data structure to hold dipole-related data and staging arrays
973 //! Dipole staging for fast summing over MPI
974 gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
975 //! Dipole staging for states A and B (index 0 and 1 resp.)
976 gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
980 static void reduceAndUpdateMuTot(DipoleData* dipoleData,
982 const bool haveFreeEnergy,
983 gmx::ArrayRef<const real> lambda,
985 const DDBalanceRegionHandler& ddBalanceRegionHandler)
989 gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
990 ddBalanceRegionHandler.reopenRegionCpu();
992 for (int i = 0; i < 2; i++)
994 for (int j = 0; j < DIM; j++)
996 dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
1000 if (!haveFreeEnergy)
1002 copy_rvec(dipoleData->muStateAB[0], muTotal);
1006 for (int j = 0; j < DIM; j++)
1008 muTotal[j] = (1.0 - lambda[efptCOUL]) * dipoleData->muStateAB[0][j]
1009 + lambda[efptCOUL] * dipoleData->muStateAB[1][j];
1014 void do_force(FILE* fplog,
1015 const t_commrec* cr,
1016 const gmx_multisim_t* ms,
1017 const t_inputrec* inputrec,
1019 gmx_enfrot* enforcedRotation,
1020 gmx::ImdSession* imdSession,
1024 gmx_wallcycle_t wcycle,
1025 const gmx_localtop_t* top,
1027 gmx::ArrayRefWithPadding<gmx::RVec> x,
1029 gmx::ForceBuffersView* forceView,
1031 const t_mdatoms* mdatoms,
1032 gmx_enerdata_t* enerd,
1033 gmx::ArrayRef<const real> lambda,
1035 gmx::MdrunScheduleWorkload* runScheduleWork,
1036 gmx::VirtualSitesHandler* vsite,
1041 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1043 auto force = forceView->forceWithPadding();
1044 GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1045 "The size of the force buffer should be at least the number of atoms to compute "
1048 nonbonded_verlet_t* nbv = fr->nbv.get();
1049 interaction_const_t* ic = fr->ic;
1050 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1052 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1055 runScheduleWork->stepWork = setupStepWorkload(legacyFlags, fr->bNonbonded, simulationWork,
1056 thisRankHasDuty(cr, DUTY_PME));
1057 const StepWorkload& stepWork = runScheduleWork->stepWork;
1060 const bool useGpuPmeOnThisRank = simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME);
1062 /* At a search step we need to start the first balancing region
1063 * somewhere early inside the step after communication during domain
1064 * decomposition (and not during the previous step as usual).
1066 if (stepWork.doNeighborSearch)
1068 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1071 clear_mat(vir_force);
1073 if (fr->pbcType != PbcType::No)
1075 /* Compute shift vectors every step,
1076 * because of pressure coupling or box deformation!
1078 if (stepWork.haveDynamicBox && stepWork.stateChanged)
1080 calc_shifts(box, fr->shift_vec);
1083 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1084 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1087 put_atoms_in_box_omp(fr->pbcType, box, x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1088 gmx_omp_nthreads_get(emntDefault));
1089 inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1093 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1095 const bool pmeSendCoordinatesFromGpu =
1096 GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1097 const bool reinitGpuPmePpComms =
1098 GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1100 const auto localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1101 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1102 AtomLocality::Local, simulationWork, stepWork)
1105 // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1106 // Otherwise the send will occur after H2D coordinate transfer.
1107 if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu)
1109 /* Send particle coordinates to the pme nodes */
1110 if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1112 GMX_RELEASE_ASSERT(false,
1113 "GPU update and separate PME ranks are only supported with GPU "
1114 "direct communication!");
1115 // TODO: when this code-path becomes supported add:
1116 // stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1119 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1120 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1121 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1122 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1125 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1126 // The local coordinates can be copied right away.
1127 // NOTE: Consider moving this copy to right after they are updated and constrained,
1128 // if the later is not offloaded.
1129 if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1131 if (stepWork.doNeighborSearch)
1133 // TODO refactor this to do_md, after partitioning.
1134 stateGpu->reinit(mdatoms->homenr,
1135 cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1136 if (useGpuPmeOnThisRank)
1138 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1139 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1142 // We need to copy coordinates when:
1143 // 1. Update is not offloaded
1144 // 2. The buffers were reinitialized on search step
1145 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1147 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1148 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1152 // TODO Update this comment when introducing SimulationWorkload
1154 // The conditions for gpuHaloExchange e.g. using GPU buffer
1155 // operations were checked before construction, so here we can
1156 // just use it and assert upon any conditions.
1157 const bool ddUsesGpuDirectCommunication =
1158 ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange.empty()));
1159 GMX_ASSERT(!ddUsesGpuDirectCommunication || stepWork.useGpuXBufferOps,
1160 "Must use coordinate buffer ops with GPU halo exchange");
1161 const bool useGpuForcesHaloExchange = ddUsesGpuDirectCommunication && stepWork.useGpuFBufferOps;
1163 // Copy coordinate from the GPU if update is on the GPU and there
1164 // are forces to be computed on the CPU, or for the computation of
1165 // virial, or if host-side data will be transferred from this task
1166 // to a remote task for halo exchange or PME-PP communication. At
1167 // search steps the current coordinates are already on the host,
1168 // hence copy is not needed.
1169 const bool haveHostPmePpComms =
1170 !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1171 const bool haveHostHaloExchangeComms = havePPDomainDecomposition(cr) && !ddUsesGpuDirectCommunication;
1173 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1174 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1175 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1176 || haveHostPmePpComms || haveHostHaloExchangeComms))
1178 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1179 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1180 haveCopiedXFromGpu = true;
1183 // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1184 // Otherwise the send will occur before the H2D coordinate transfer.
1185 if (!thisRankHasDuty(cr, DUTY_PME) && pmeSendCoordinatesFromGpu)
1187 /* Send particle coordinates to the pme nodes */
1188 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1189 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1190 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1191 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1194 if (useGpuPmeOnThisRank)
1196 launchPmeGpuSpread(fr->pmedata, box, stepWork, localXReadyOnDevice, lambda[efptCOUL], wcycle);
1199 /* do gridding for pair search */
1200 if (stepWork.doNeighborSearch)
1202 if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1204 fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1208 // - vzero is constant, do we need to pass it?
1209 // - box_diag should be passed directly to nbnxn_put_on_grid
1215 box_diag[XX] = box[XX][XX];
1216 box_diag[YY] = box[YY][YY];
1217 box_diag[ZZ] = box[ZZ][ZZ];
1219 wallcycle_start(wcycle, ewcNS);
1220 if (!DOMAINDECOMP(cr))
1222 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1223 nbnxn_put_on_grid(nbv, box, 0, vzero, box_diag, nullptr, { 0, mdatoms->homenr }, -1,
1224 fr->cginfo, x.unpaddedArrayRef(), 0, nullptr);
1225 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1229 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1230 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1231 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1234 nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1235 gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr), fr->cginfo);
1237 wallcycle_stop(wcycle, ewcNS);
1239 /* initialize the GPU nbnxm atom data and bonded data structures */
1240 if (simulationWork.useGpuNonbonded)
1242 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1244 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1245 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1246 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1250 /* Now we put all atoms on the grid, we can assign bonded
1251 * interactions to the GPU, where the grid order is
1252 * needed. Also the xq, f and fshift device buffers have
1253 * been reallocated if needed, so the bonded code can
1254 * learn about them. */
1255 // TODO the xq, f, and fshift buffers are now shared
1256 // resources, so they should be maintained by a
1257 // higher-level object than the nb module.
1258 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(
1259 nbv->getGridIndices(), top->idef, Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1260 Nbnxm::gpu_get_f(nbv->gpu_nbv), Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1262 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1265 // Need to run after the GPU-offload bonded interaction lists
1266 // are set up to be able to determine whether there is bonded work.
1267 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1268 *inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1270 wallcycle_start_nocount(wcycle, ewcNS);
1271 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1272 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1273 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1275 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1277 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1278 wallcycle_stop(wcycle, ewcNS);
1280 if (stepWork.useGpuXBufferOps)
1282 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1284 // For force buffer ops, we use the below conditon rather than
1285 // useGpuFBufferOps to ensure that init is performed even if this
1286 // NS step is also a virial step (on which f buf ops are deactivated).
1287 if (GMX_GPU_CUDA && simulationWork.useGpuBufferOps && simulationWork.useGpuNonbonded)
1289 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1290 nbv->atomdata_init_add_nbat_f_to_f_gpu(stateGpu->fReducedOnDevice());
1293 else if (!EI_TPI(inputrec->eI))
1295 if (stepWork.useGpuXBufferOps)
1297 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1298 nbv->convertCoordinatesGpu(AtomLocality::Local, false, stateGpu->getCoordinates(),
1299 localXReadyOnDevice);
1303 if (simulationWork.useGpuUpdate)
1305 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1306 GMX_ASSERT(haveCopiedXFromGpu,
1307 "a wait should only be triggered if copy has been scheduled");
1308 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1310 nbv->convertCoordinates(AtomLocality::Local, false, x.unpaddedArrayRef());
1314 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1316 if (simulationWork.useGpuNonbonded)
1318 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1320 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1322 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1323 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1324 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1326 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1328 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1329 // with X buffer ops offloaded to the GPU on all but the search steps
1331 // bonded work not split into separate local and non-local, so with DD
1332 // we can only launch the kernel after non-local coordinates have been received.
1333 if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1335 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1336 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1337 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1340 /* launch local nonbonded work on GPU */
1341 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1342 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1343 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1344 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1347 if (useGpuPmeOnThisRank)
1349 // In PME GPU and mixed mode we launch FFT / gather after the
1350 // X copy/transform to allow overlap as well as after the GPU NB
1351 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1352 // the nonbonded kernel.
1353 launchPmeGpuFftAndGather(fr->pmedata, lambda[efptCOUL], wcycle, stepWork);
1356 /* Communicate coordinates and sum dipole if necessary +
1357 do non-local pair search */
1358 if (havePPDomainDecomposition(cr))
1360 if (stepWork.doNeighborSearch)
1362 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1363 wallcycle_start_nocount(wcycle, ewcNS);
1364 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1365 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1366 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1368 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1369 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1370 wallcycle_stop(wcycle, ewcNS);
1371 // TODO refactor this GPU halo exchange re-initialisation
1372 // to location in do_md where GPU halo exchange is
1373 // constructed at partitioning, after above stateGpu
1374 // re-initialization has similarly been refactored
1375 if (ddUsesGpuDirectCommunication)
1377 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1382 if (ddUsesGpuDirectCommunication)
1384 // The following must be called after local setCoordinates (which records an event
1385 // when the coordinate data has been copied to the device).
1386 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1388 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1390 // non-local part of coordinate buffer must be copied back to host for CPU work
1391 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1396 // Note: GPU update + DD without direct communication is not supported,
1397 // a waitCoordinatesReadyOnHost() should be issued if it will be.
1398 GMX_ASSERT(!simulationWork.useGpuUpdate,
1399 "GPU update is not supported with CPU halo exchange");
1400 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1403 if (stepWork.useGpuXBufferOps)
1405 if (!useGpuPmeOnThisRank && !ddUsesGpuDirectCommunication)
1407 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1409 nbv->convertCoordinatesGpu(AtomLocality::NonLocal, false, stateGpu->getCoordinates(),
1410 stateGpu->getCoordinatesReadyOnDeviceEvent(
1411 AtomLocality::NonLocal, simulationWork, stepWork));
1415 nbv->convertCoordinates(AtomLocality::NonLocal, false, x.unpaddedArrayRef());
1419 if (simulationWork.useGpuNonbonded)
1421 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1423 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1425 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1426 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1427 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1430 if (domainWork.haveGpuBondedWork)
1432 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1433 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1434 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1437 /* launch non-local nonbonded tasks on GPU */
1438 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1439 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1441 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1443 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1447 if (simulationWork.useGpuNonbonded)
1449 /* launch D2H copy-back F */
1450 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1451 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1453 if (havePPDomainDecomposition(cr))
1455 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1457 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1458 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1460 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1462 fr->gpuBonded->launchEnergyTransfer();
1464 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1467 gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1468 if (fr->wholeMoleculeTransform)
1470 xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1473 DipoleData dipoleData;
1475 if (simulationWork.computeMuTot)
1477 const int start = 0;
1479 /* Calculate total (local) dipole moment in a temporary common array.
1480 * This makes it possible to sum them over nodes faster.
1482 gmx::ArrayRef<const gmx::RVec> xRef =
1483 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1484 calc_mu(start, mdatoms->homenr, xRef, mdatoms->chargeA, mdatoms->chargeB,
1485 mdatoms->nChargePerturbed, dipoleData.muStaging[0], dipoleData.muStaging[1]);
1487 reduceAndUpdateMuTot(&dipoleData, cr, (fr->efep != efepNO), lambda, muTotal, ddBalanceRegionHandler);
1490 /* Reset energies */
1491 reset_enerdata(enerd);
1493 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1495 wallcycle_start(wcycle, ewcPPDURINGPME);
1496 dd_force_flop_start(cr->dd, nrnb);
1499 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1500 // this wait ensures that the D2H transfer is complete.
1501 if ((simulationWork.useGpuUpdate)
1502 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1504 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1509 wallcycle_start(wcycle, ewcROT);
1510 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step,
1511 stepWork.doNeighborSearch);
1512 wallcycle_stop(wcycle, ewcROT);
1515 /* Start the force cycle counter.
1516 * Note that a different counter is used for dynamic load balancing.
1518 wallcycle_start(wcycle, ewcFORCE);
1520 // Set up and clear force outputs.
1521 // We use std::move to keep the compiler happy, it has no effect.
1522 ForceOutputs forceOut = setupForceOutputs(fr->forceHelperBuffers.get(), pull_work, *inputrec,
1523 std::move(force), stepWork, wcycle);
1525 /* We calculate the non-bonded forces, when done on the CPU, here.
1526 * We do this before calling do_force_lowlevel, because in that
1527 * function, the listed forces are calculated before PME, which
1528 * does communication. With this order, non-bonded and listed
1529 * force calculation imbalance can be balanced out by the domain
1530 * decomposition load balancing.
1533 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1535 if (!useOrEmulateGpuNb)
1537 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1540 if (fr->efep != efepNO)
1542 /* Calculate the local and non-local free energy interactions here.
1543 * Happens here on the CPU both with and without GPU.
1545 nbv->dispatchFreeEnergyKernel(InteractionLocality::Local, fr,
1546 as_rvec_array(x.unpaddedArrayRef().data()),
1547 &forceOut.forceWithShiftForces(), *mdatoms, inputrec->fepvals,
1548 lambda, enerd, stepWork, nrnb);
1550 if (havePPDomainDecomposition(cr))
1552 nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal, fr,
1553 as_rvec_array(x.unpaddedArrayRef().data()),
1554 &forceOut.forceWithShiftForces(), *mdatoms,
1555 inputrec->fepvals, lambda, enerd, stepWork, nrnb);
1559 if (!useOrEmulateGpuNb)
1561 if (havePPDomainDecomposition(cr))
1563 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1567 if (stepWork.computeForces)
1569 /* Add all the non-bonded force to the normal force array.
1570 * This can be split into a local and a non-local part when overlapping
1571 * communication with calculation with domain decomposition.
1573 wallcycle_stop(wcycle, ewcFORCE);
1574 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All, forceOut.forceWithShiftForces().force());
1575 wallcycle_start_nocount(wcycle, ewcFORCE);
1578 /* If there are multiple fshift output buffers we need to reduce them */
1579 if (stepWork.computeVirial)
1581 /* This is not in a subcounter because it takes a
1582 negligible and constant-sized amount of time */
1583 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat,
1584 forceOut.forceWithShiftForces().shiftForces());
1588 // TODO Force flags should include haveFreeEnergyWork for this domain
1589 if (ddUsesGpuDirectCommunication && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1591 /* Wait for non-local coordinate data to be copied from device */
1592 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1594 /* Compute the bonded and non-bonded energies and optionally forces */
1595 do_force_lowlevel(fr, inputrec, cr, ms, nrnb, wcycle, mdatoms, x, xWholeMolecules, hist,
1596 &forceOut, enerd, box, lambda.data(), as_rvec_array(dipoleData.muStateAB),
1597 stepWork, ddBalanceRegionHandler);
1599 wallcycle_stop(wcycle, ewcFORCE);
1601 // VdW dispersion correction, only computed on master rank to avoid double counting
1602 if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1604 // Calculate long range corrections to pressure and energy
1605 const DispersionCorrection::Correction correction =
1606 fr->dispersionCorrection->calculate(box, lambda[efptVDW]);
1608 if (stepWork.computeEnergy)
1610 enerd->term[F_DISPCORR] = correction.energy;
1611 enerd->term[F_DVDL_VDW] += correction.dvdl;
1612 enerd->dvdl_lin[efptVDW] += correction.dvdl;
1614 if (stepWork.computeVirial)
1616 correction.correctVirial(vir_force);
1617 enerd->term[F_PDISPCORR] = correction.pressure;
1621 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation, imdSession, pull_work, step, t,
1622 wcycle, fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda,
1623 stepWork, &forceOut.forceWithVirial(), enerd, ed, stepWork.doNeighborSearch);
1626 // Will store the amount of cycles spent waiting for the GPU that
1627 // will be later used in the DLB accounting.
1628 float cycles_wait_gpu = 0;
1629 if (useOrEmulateGpuNb)
1631 auto& forceWithShiftForces = forceOut.forceWithShiftForces();
1633 /* wait for non-local forces (or calculate in emulation mode) */
1634 if (havePPDomainDecomposition(cr))
1636 if (simulationWork.useGpuNonbonded)
1638 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1639 nbv->gpu_nbv, stepWork, AtomLocality::NonLocal, enerd->grpp.ener[egLJSR].data(),
1640 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(), wcycle);
1644 wallcycle_start_nocount(wcycle, ewcFORCE);
1645 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes,
1646 step, nrnb, wcycle);
1647 wallcycle_stop(wcycle, ewcFORCE);
1650 if (stepWork.useGpuFBufferOps)
1652 gmx::FixedCapacityVector<GpuEventSynchronizer*, 1> dependencyList;
1654 // TODO: move this into DomainLifetimeWorkload, including the second part of the
1655 // condition The bonded and free energy CPU tasks can have non-local force
1656 // contributions which are a dependency for the GPU force reduction.
1657 bool haveNonLocalForceContribInCpuBuffer =
1658 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1660 if (haveNonLocalForceContribInCpuBuffer)
1662 stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(),
1663 AtomLocality::NonLocal);
1664 dependencyList.push_back(stateGpu->getForcesReadyOnDeviceEvent(
1665 AtomLocality::NonLocal, stepWork.useGpuFBufferOps));
1668 nbv->atomdata_add_nbat_f_to_f_gpu(AtomLocality::NonLocal, stateGpu->getForces(),
1669 pme_gpu_get_device_f(fr->pmedata), dependencyList,
1670 false, haveNonLocalForceContribInCpuBuffer);
1671 if (!useGpuForcesHaloExchange)
1673 // copy from GPU input for dd_move_f()
1674 stateGpu->copyForcesFromGpu(forceOut.forceWithShiftForces().force(),
1675 AtomLocality::NonLocal);
1680 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
1684 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
1686 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
1691 if (havePPDomainDecomposition(cr))
1693 /* We are done with the CPU compute.
1694 * We will now communicate the non-local forces.
1695 * If we use a GPU this will overlap with GPU work, so in that case
1696 * we do not close the DD force balancing region here.
1698 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1700 if (stepWork.computeForces)
1703 if (useGpuForcesHaloExchange)
1705 if (domainWork.haveCpuLocalForceWork)
1707 stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(), AtomLocality::Local);
1709 communicateGpuHaloForces(*cr, domainWork.haveCpuLocalForceWork);
1713 if (stepWork.useGpuFBufferOps)
1715 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
1717 dd_move_f(cr->dd, &forceOut.forceWithShiftForces(), wcycle);
1722 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1723 // an alternating wait/reduction scheme.
1724 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
1725 && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
1726 if (alternateGpuWait)
1728 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &forceOut, enerd, lambda[efptCOUL],
1732 if (!alternateGpuWait && useGpuPmeOnThisRank)
1734 pme_gpu_wait_and_reduce(fr->pmedata, stepWork, wcycle, &forceOut.forceWithVirial(), enerd,
1738 /* Wait for local GPU NB outputs on the non-alternating wait path */
1739 if (!alternateGpuWait && simulationWork.useGpuNonbonded)
1741 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1742 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1743 * but even with a step of 0.1 ms the difference is less than 1%
1746 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
1747 const float waitCycles = Nbnxm::gpu_wait_finish_task(
1748 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
1749 enerd->grpp.ener[egCOULSR].data(), forceOut.forceWithShiftForces().shiftForces(), wcycle);
1751 if (ddBalanceRegionHandler.useBalancingRegion())
1753 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1754 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
1756 /* We measured few cycles, it could be that the kernel
1757 * and transfer finished earlier and there was no actual
1758 * wait time, only API call overhead.
1759 * Then the actual time could be anywhere between 0 and
1760 * cycles_wait_est. We will use half of cycles_wait_est.
1762 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1764 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1768 if (fr->nbv->emulateGpu())
1770 // NOTE: emulation kernel is not included in the balancing region,
1771 // but emulation mode does not target performance anyway
1772 wallcycle_start_nocount(wcycle, ewcFORCE);
1773 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local,
1774 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes, step, nrnb, wcycle);
1775 wallcycle_stop(wcycle, ewcFORCE);
1778 // If on GPU PME-PP comms or GPU update path, receive forces from PME before GPU buffer ops
1779 // TODO refactor this and unify with below default-path call to the same function
1780 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME)
1781 && (simulationWork.useGpuPmePpCommunication || simulationWork.useGpuUpdate))
1783 /* In case of node-splitting, the PP nodes receive the long-range
1784 * forces, virial and energy from the PME nodes here.
1786 pme_receive_force_ener(fr, cr, &forceOut.forceWithVirial(), enerd,
1787 simulationWork.useGpuPmePpCommunication,
1788 stepWork.useGpuPmeFReduction, wcycle);
1792 /* Do the nonbonded GPU (or emulation) force buffer reduction
1793 * on the non-alternating path. */
1794 if (useOrEmulateGpuNb && !alternateGpuWait)
1796 // TODO simplify the below conditionals. Pass buffer and sync pointers at init stage rather than here. Unify getter fns for sameGPU/otherGPU cases.
1798 stepWork.useGpuPmeFReduction
1799 ? (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1800 : // PME force buffer on same GPU
1801 fr->pmePpCommGpu->getGpuForceStagingPtr()) // buffer received from other GPU
1802 : nullptr; // PME reduction not active on GPU
1804 GpuEventSynchronizer* const pmeSynchronizer =
1805 stepWork.useGpuPmeFReduction
1806 ? (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1807 : // PME force buffer on same GPU
1808 static_cast<GpuEventSynchronizer*>(
1809 fr->pmePpCommGpu->getForcesReadySynchronizer())) // buffer received from other GPU
1810 : nullptr; // PME reduction not active on GPU
1812 gmx::FixedCapacityVector<GpuEventSynchronizer*, 3> dependencyList;
1814 if (stepWork.useGpuPmeFReduction)
1816 dependencyList.push_back(pmeSynchronizer);
1819 gmx::ArrayRef<gmx::RVec> forceWithShift = forceOut.forceWithShiftForces().force();
1821 if (stepWork.useGpuFBufferOps)
1823 // Flag to specify whether the CPU force buffer has contributions to
1824 // local atoms. This depends on whether there are CPU-based force tasks
1825 // or when DD is active the halo exchange has resulted in contributions
1826 // from the non-local part.
1827 const bool haveLocalForceContribInCpuBuffer =
1828 (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
1830 // TODO: move these steps as early as possible:
1831 // - CPU f H2D should be as soon as all CPU-side forces are done
1832 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
1833 // before the next CPU task that consumes the forces: vsite spread or update)
1834 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
1835 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
1836 // These should be unified.
1837 if (haveLocalForceContribInCpuBuffer && !useGpuForcesHaloExchange)
1839 // Note: AtomLocality::All is used for the non-DD case because, as in this
1840 // case copyForcesToGpu() uses a separate stream, it allows overlap of
1841 // CPU force H2D with GPU force tasks on all streams including those in the
1842 // local stream which would otherwise be implicit dependencies for the
1843 // transfer and would not overlap.
1844 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
1846 stateGpu->copyForcesToGpu(forceWithShift, locality);
1847 dependencyList.push_back(
1848 stateGpu->getForcesReadyOnDeviceEvent(locality, stepWork.useGpuFBufferOps));
1850 if (useGpuForcesHaloExchange)
1852 dependencyList.push_back(cr->dd->gpuHaloExchange[0]->getForcesReadyOnDeviceEvent());
1854 nbv->atomdata_add_nbat_f_to_f_gpu(AtomLocality::Local, stateGpu->getForces(), pmeForcePtr,
1855 dependencyList, stepWork.useGpuPmeFReduction,
1856 haveLocalForceContribInCpuBuffer);
1857 // Copy forces to host if they are needed for update or if virtual sites are enabled.
1858 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
1859 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1860 // copy call done in sim_utils(...) for the output.
1861 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
1862 // they should not be copied in do_md(...) for the output.
1863 if (!simulationWork.useGpuUpdate || vsite)
1865 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
1866 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1871 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
1875 launchGpuEndOfStepTasks(nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork,
1876 useGpuPmeOnThisRank, step, wcycle);
1878 if (DOMAINDECOMP(cr))
1880 dd_force_flop_stop(cr->dd, nrnb);
1883 if (stepWork.computeForces)
1885 postProcessForceWithShiftForces(nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOut,
1886 vir_force, *mdatoms, *fr, vsite, stepWork);
1889 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
1890 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
1891 && !simulationWork.useGpuUpdate)
1893 /* In case of node-splitting, the PP nodes receive the long-range
1894 * forces, virial and energy from the PME nodes here.
1896 pme_receive_force_ener(fr, cr, &forceOut.forceWithVirial(), enerd,
1897 simulationWork.useGpuPmePpCommunication, false, wcycle);
1900 if (stepWork.computeForces)
1902 postProcessForces(cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOut, vir_force,
1903 mdatoms, fr, vsite, stepWork);
1906 if (stepWork.computeEnergy)
1908 /* Compute the final potential energy terms */
1909 accumulatePotentialEnergies(enerd, lambda, inputrec->fepvals);
1911 if (!EI_TPI(inputrec->eI))
1913 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1917 /* In case we don't have constraints and are using GPUs, the next balancing
1918 * region starts here.
1919 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1920 * virial calculation and COM pulling, is not thus not included in
1921 * the balance timing, which is ok as most tasks do communication.
1923 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);