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/wall.h"
85 #include "gromacs/mdlib/wholemoleculetransform.h"
86 #include "gromacs/mdtypes/commrec.h"
87 #include "gromacs/mdtypes/enerdata.h"
88 #include "gromacs/mdtypes/forcebuffers.h"
89 #include "gromacs/mdtypes/forceoutput.h"
90 #include "gromacs/mdtypes/forcerec.h"
91 #include "gromacs/mdtypes/iforceprovider.h"
92 #include "gromacs/mdtypes/inputrec.h"
93 #include "gromacs/mdtypes/md_enums.h"
94 #include "gromacs/mdtypes/mdatom.h"
95 #include "gromacs/mdtypes/simulation_workload.h"
96 #include "gromacs/mdtypes/state.h"
97 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
98 #include "gromacs/nbnxm/gpu_data_mgmt.h"
99 #include "gromacs/nbnxm/nbnxm.h"
100 #include "gromacs/nbnxm/nbnxm_gpu.h"
101 #include "gromacs/pbcutil/ishift.h"
102 #include "gromacs/pbcutil/pbc.h"
103 #include "gromacs/pulling/pull.h"
104 #include "gromacs/pulling/pull_rotation.h"
105 #include "gromacs/timing/cyclecounter.h"
106 #include "gromacs/timing/gpu_timing.h"
107 #include "gromacs/timing/wallcycle.h"
108 #include "gromacs/timing/wallcyclereporting.h"
109 #include "gromacs/timing/walltime_accounting.h"
110 #include "gromacs/topology/topology.h"
111 #include "gromacs/utility/arrayref.h"
112 #include "gromacs/utility/basedefinitions.h"
113 #include "gromacs/utility/cstringutil.h"
114 #include "gromacs/utility/exceptions.h"
115 #include "gromacs/utility/fatalerror.h"
116 #include "gromacs/utility/fixedcapacityvector.h"
117 #include "gromacs/utility/gmxassert.h"
118 #include "gromacs/utility/gmxmpi.h"
119 #include "gromacs/utility/logger.h"
120 #include "gromacs/utility/smalloc.h"
121 #include "gromacs/utility/strconvert.h"
122 #include "gromacs/utility/sysinfo.h"
125 using gmx::AtomLocality;
126 using gmx::DomainLifetimeWorkload;
127 using gmx::ForceOutputs;
128 using gmx::ForceWithShiftForces;
129 using gmx::InteractionLocality;
131 using gmx::SimulationWorkload;
132 using gmx::StepWorkload;
134 // TODO: this environment variable allows us to verify before release
135 // that on less common architectures the total cost of polling is not larger than
136 // a blocking wait (so polling does not introduce overhead when the static
137 // PME-first ordering would suffice).
138 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
140 static void sum_forces(ArrayRef<RVec> f, ArrayRef<const RVec> forceToAdd)
142 GMX_ASSERT(f.size() >= forceToAdd.size(), "Accumulation buffer should be sufficiently large");
143 const int end = forceToAdd.size();
145 int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
146 #pragma omp parallel for num_threads(nt) schedule(static)
147 for (int i = 0; i < end; i++)
149 rvec_inc(f[i], forceToAdd[i]);
153 static void calc_virial(int start,
156 const gmx::ForceWithShiftForces& forceWithShiftForces,
160 const t_forcerec* fr,
163 /* The short-range virial from surrounding boxes */
164 const rvec* fshift = as_rvec_array(forceWithShiftForces.shiftForces().data());
165 calc_vir(SHIFTS, fr->shift_vec, fshift, vir_part, pbcType == PbcType::Screw, box);
166 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
168 /* Calculate partial virial, for local atoms only, based on short range.
169 * Total virial is computed in global_stat, called from do_md
171 const rvec* f = as_rvec_array(forceWithShiftForces.force().data());
172 f_calc_vir(start, start + homenr, x, f, vir_part, box);
173 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
177 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
181 static void pull_potential_wrapper(const t_commrec* cr,
182 const t_inputrec* ir,
184 gmx::ArrayRef<const gmx::RVec> x,
185 gmx::ForceWithVirial* force,
186 const t_mdatoms* mdatoms,
187 gmx_enerdata_t* enerd,
191 gmx_wallcycle_t wcycle)
196 /* Calculate the center of mass forces, this requires communication,
197 * which is why pull_potential is called close to other communication.
199 wallcycle_start(wcycle, ewcPULLPOT);
200 set_pbc(&pbc, ir->pbcType, box);
202 enerd->term[F_COM_PULL] +=
203 pull_potential(pull_work, mdatoms->massT, &pbc, cr, t, lambda[efptRESTRAINT],
204 as_rvec_array(x.data()), force, &dvdl);
205 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
206 wallcycle_stop(wcycle, ewcPULLPOT);
209 static void pme_receive_force_ener(t_forcerec* fr,
211 gmx::ForceWithVirial* forceWithVirial,
212 gmx_enerdata_t* enerd,
213 bool useGpuPmePpComms,
214 bool receivePmeForceToGpu,
215 gmx_wallcycle_t wcycle)
217 real e_q, e_lj, dvdl_q, dvdl_lj;
218 float cycles_ppdpme, cycles_seppme;
220 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
221 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
223 /* In case of node-splitting, the PP nodes receive the long-range
224 * forces, virial and energy from the PME nodes here.
226 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
229 gmx_pme_receive_f(fr->pmePpCommGpu.get(), cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
230 useGpuPmePpComms, receivePmeForceToGpu, &cycles_seppme);
231 enerd->term[F_COUL_RECIP] += e_q;
232 enerd->term[F_LJ_RECIP] += e_lj;
233 enerd->dvdl_lin[efptCOUL] += dvdl_q;
234 enerd->dvdl_lin[efptVDW] += dvdl_lj;
238 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
240 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
243 static void print_large_forces(FILE* fp,
248 ArrayRef<const RVec> x,
249 ArrayRef<const RVec> f)
251 real force2Tolerance = gmx::square(forceTolerance);
252 gmx::index numNonFinite = 0;
253 for (int i = 0; i < md->homenr; i++)
255 real force2 = norm2(f[i]);
256 bool nonFinite = !std::isfinite(force2);
257 if (force2 >= force2Tolerance || nonFinite)
259 fprintf(fp, "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n", step,
260 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
267 if (numNonFinite > 0)
269 /* Note that with MPI this fatal call on one rank might interrupt
270 * the printing on other ranks. But we can only avoid that with
271 * an expensive MPI barrier that we would need at each step.
273 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
277 //! When necessary, spreads forces on vsites and computes the virial for \p forceOutputs->forceWithShiftForces()
278 static void postProcessForceWithShiftForces(t_nrnb* nrnb,
279 gmx_wallcycle_t wcycle,
281 ArrayRef<const RVec> x,
282 ForceOutputs* forceOutputs,
284 const t_mdatoms& mdatoms,
285 const t_forcerec& fr,
286 gmx::VirtualSitesHandler* vsite,
287 const StepWorkload& stepWork)
289 ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
291 /* If we have NoVirSum forces, but we do not calculate the virial,
292 * we later sum the forceWithShiftForces buffer together with
293 * the noVirSum buffer and spread the combined vsite forces at once.
295 if (vsite && (!forceOutputs->haveForceWithVirial() || stepWork.computeVirial))
297 using VirialHandling = gmx::VirtualSitesHandler::VirialHandling;
299 auto f = forceWithShiftForces.force();
300 auto fshift = forceWithShiftForces.shiftForces();
301 const VirialHandling virialHandling =
302 (stepWork.computeVirial ? VirialHandling::Pbc : VirialHandling::None);
303 vsite->spreadForces(x, f, virialHandling, fshift, nullptr, nrnb, box, wcycle);
304 forceWithShiftForces.haveSpreadVsiteForces() = true;
307 if (stepWork.computeVirial)
309 /* Calculation of the virial must be done after vsites! */
310 calc_virial(0, mdatoms.homenr, as_rvec_array(x.data()), forceWithShiftForces, vir_force,
311 box, nrnb, &fr, fr.pbcType);
315 //! Spread, compute virial for and sum forces, when necessary
316 static void postProcessForces(const t_commrec* cr,
319 gmx_wallcycle_t wcycle,
321 ArrayRef<const RVec> x,
322 ForceOutputs* forceOutputs,
324 const t_mdatoms* mdatoms,
325 const t_forcerec* fr,
326 gmx::VirtualSitesHandler* vsite,
327 const StepWorkload& stepWork)
329 // Extract the final output force buffer, which is also the buffer for forces with shift forces
330 ArrayRef<RVec> f = forceOutputs->forceWithShiftForces().force();
332 if (forceOutputs->haveForceWithVirial())
334 auto& forceWithVirial = forceOutputs->forceWithVirial();
338 /* Spread the mesh force on virtual sites to the other particles...
339 * This is parallellized. MPI communication is performed
340 * if the constructing atoms aren't local.
342 GMX_ASSERT(!stepWork.computeVirial || f.data() != forceWithVirial.force_.data(),
343 "We need separate force buffers for shift and virial forces when "
344 "computing the virial");
345 GMX_ASSERT(!stepWork.computeVirial
346 || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
347 "We should spread the force with shift forces separately when computing "
349 const gmx::VirtualSitesHandler::VirialHandling virialHandling =
350 (stepWork.computeVirial ? gmx::VirtualSitesHandler::VirialHandling::NonLinear
351 : gmx::VirtualSitesHandler::VirialHandling::None);
352 matrix virial = { { 0 } };
353 vsite->spreadForces(x, forceWithVirial.force_, virialHandling, {}, virial, nrnb, box, wcycle);
354 forceWithVirial.addVirialContribution(virial);
357 if (stepWork.computeVirial)
359 /* Now add the forces, this is local */
360 sum_forces(f, forceWithVirial.force_);
362 /* Add the direct virial contributions */
364 forceWithVirial.computeVirial_,
365 "forceWithVirial should request virial computation when we request the virial");
366 m_add(vir_force, forceWithVirial.getVirial(), vir_force);
370 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
376 GMX_ASSERT(vsite == nullptr || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
377 "We should have spread the vsite forces (earlier)");
380 if (fr->print_force >= 0)
382 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
386 static void do_nb_verlet(t_forcerec* fr,
387 const interaction_const_t* ic,
388 gmx_enerdata_t* enerd,
389 const StepWorkload& stepWork,
390 const InteractionLocality ilocality,
394 gmx_wallcycle_t wcycle)
396 if (!stepWork.computeNonbondedForces)
398 /* skip non-bonded calculation */
402 nonbonded_verlet_t* nbv = fr->nbv.get();
404 /* GPU kernel launch overhead is already timed separately */
405 if (fr->cutoff_scheme != ecutsVERLET)
407 gmx_incons("Invalid cut-off scheme passed!");
412 /* When dynamic pair-list pruning is requested, we need to prune
413 * at nstlistPrune steps.
415 if (nbv->isDynamicPruningStepCpu(step))
417 /* Prune the pair-list beyond fr->ic->rlistPrune using
418 * the current coordinates of the atoms.
420 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
421 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
422 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
426 nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
429 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
431 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, v.ssize());
433 /* Note that we would like to avoid this conditional by putting it
434 * into the omp pragma instead, but then we still take the full
435 * omp parallel for overhead (at least with gcc5).
437 if (!useOpenmpThreading || nth == 1)
446 #pragma omp parallel for num_threads(nth) schedule(static)
447 for (gmx::index i = 0; i < v.ssize(); i++)
454 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
456 * \param groupOptions Group options, containing T-coupling options
458 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
460 real nrdfCoupled = 0;
461 real nrdfUncoupled = 0;
462 real kineticEnergy = 0;
463 for (int g = 0; g < groupOptions.ngtc; g++)
465 if (groupOptions.tau_t[g] >= 0)
467 nrdfCoupled += groupOptions.nrdf[g];
468 kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * BOLTZ;
472 nrdfUncoupled += groupOptions.nrdf[g];
476 /* This conditional with > also catches nrdf=0 */
477 if (nrdfCoupled > nrdfUncoupled)
479 return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
487 /*! \brief This routine checks that the potential energy is finite.
489 * Always checks that the potential energy is finite. If step equals
490 * inputrec.init_step also checks that the magnitude of the potential energy
491 * is reasonable. Terminates with a fatal error when a check fails.
492 * Note that passing this check does not guarantee finite forces,
493 * since those use slightly different arithmetics. But in most cases
494 * there is just a narrow coordinate range where forces are not finite
495 * and energies are finite.
497 * \param[in] step The step number, used for checking and printing
498 * \param[in] enerd The energy data; the non-bonded group energies need to be added to
499 * enerd.term[F_EPOT] before calling this routine \param[in] inputrec The input record
501 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
503 /* Threshold valid for comparing absolute potential energy against
504 * the kinetic energy. Normally one should not consider absolute
505 * potential energy values, but with a factor of one million
506 * we should never get false positives.
508 constexpr real c_thresholdFactor = 1e6;
510 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
511 real averageKineticEnergy = 0;
512 /* We only check for large potential energy at the initial step,
513 * because that is by far the most likely step for this too occur
514 * and because computing the average kinetic energy is not free.
515 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
516 * before they become NaN.
518 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
520 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
523 if (energyIsNotFinite
524 || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
529 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
530 "contributions to the energy are %g and %g, respectively. A %s potential energy "
531 "can be caused by overlapping interactions in bonded interactions or very large%s "
532 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
533 "configuration, incorrect interactions or parameters in the topology.",
534 step, enerd.term[F_EPOT], energyIsNotFinite ? "not finite" : "extremely high",
535 enerd.term[F_LJ], enerd.term[F_COUL_SR],
536 energyIsNotFinite ? "non-finite" : "very high", energyIsNotFinite ? " or Nan" : "");
540 /*! \brief Return true if there are special forces computed this step.
542 * The conditionals exactly correspond to those in computeSpecialForces().
544 static bool haveSpecialForces(const t_inputrec& inputrec,
545 const gmx::ForceProviders& forceProviders,
546 const pull_t* pull_work,
547 const bool computeForces,
551 return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
552 (inputrec.bPull && pull_have_potential(pull_work)) || // pull
553 inputrec.bRot || // enforced rotation
554 (ed != nullptr) || // flooding
555 (inputrec.bIMD && computeForces)); // IMD
558 /*! \brief Compute forces and/or energies for special algorithms
560 * The intention is to collect all calls to algorithms that compute
561 * forces on local atoms only and that do not contribute to the local
562 * virial sum (but add their virial contribution separately).
563 * Eventually these should likely all become ForceProviders.
564 * Within this function the intention is to have algorithms that do
565 * global communication at the end, so global barriers within the MD loop
566 * are as close together as possible.
568 * \param[in] fplog The log file
569 * \param[in] cr The communication record
570 * \param[in] inputrec The input record
571 * \param[in] awh The Awh module (nullptr if none in use).
572 * \param[in] enforcedRotation Enforced rotation module.
573 * \param[in] imdSession The IMD session
574 * \param[in] pull_work The pull work structure.
575 * \param[in] step The current MD step
576 * \param[in] t The current time
577 * \param[in,out] wcycle Wallcycle accounting struct
578 * \param[in,out] forceProviders Pointer to a list of force providers
579 * \param[in] box The unit cell
580 * \param[in] x The coordinates
581 * \param[in] mdatoms Per atom properties
582 * \param[in] lambda Array of free-energy lambda values
583 * \param[in] stepWork Step schedule flags
584 * \param[in,out] forceWithVirial Force and virial buffers
585 * \param[in,out] enerd Energy buffer
586 * \param[in,out] ed Essential dynamics pointer
587 * \param[in] didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
589 * \todo Remove didNeighborSearch, which is used incorrectly.
590 * \todo Convert all other algorithms called here to ForceProviders.
592 static void computeSpecialForces(FILE* fplog,
594 const t_inputrec* inputrec,
596 gmx_enfrot* enforcedRotation,
597 gmx::ImdSession* imdSession,
601 gmx_wallcycle_t wcycle,
602 gmx::ForceProviders* forceProviders,
604 gmx::ArrayRef<const gmx::RVec> x,
605 const t_mdatoms* mdatoms,
606 gmx::ArrayRef<const real> lambda,
607 const StepWorkload& stepWork,
608 gmx::ForceWithVirial* forceWithVirial,
609 gmx_enerdata_t* enerd,
611 bool didNeighborSearch)
613 /* NOTE: Currently all ForceProviders only provide forces.
614 * When they also provide energies, remove this conditional.
616 if (stepWork.computeForces)
618 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
619 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
621 /* Collect forces from modules */
622 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
625 if (inputrec->bPull && pull_have_potential(pull_work))
627 pull_potential_wrapper(cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work,
628 lambda.data(), t, wcycle);
632 const bool needForeignEnergyDifferences = awh->needForeignEnergyDifferences(step);
633 std::vector<double> foreignLambdaDeltaH, foreignLambdaDhDl;
634 if (needForeignEnergyDifferences)
636 enerd->foreignLambdaTerms.finalizePotentialContributions(enerd->dvdl_lin, lambda,
638 std::tie(foreignLambdaDeltaH, foreignLambdaDhDl) = enerd->foreignLambdaTerms.getTerms(cr);
641 enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(
642 inputrec->pbcType, mdatoms->massT, foreignLambdaDeltaH, foreignLambdaDhDl, box,
643 forceWithVirial, t, step, wcycle, fplog);
646 rvec* f = as_rvec_array(forceWithVirial->force_.data());
648 /* Add the forces from enforced rotation potentials (if any) */
651 wallcycle_start(wcycle, ewcROTadd);
652 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
653 wallcycle_stop(wcycle, ewcROTadd);
658 /* Note that since init_edsam() is called after the initialization
659 * of forcerec, edsam doesn't request the noVirSum force buffer.
660 * Thus if no other algorithm (e.g. PME) requires it, the forces
661 * here will contribute to the virial.
663 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
666 /* Add forces from interactive molecular dynamics (IMD), if any */
667 if (inputrec->bIMD && stepWork.computeForces)
669 imdSession->applyForces(f);
673 /*! \brief Launch the prepare_step and spread stages of PME GPU.
675 * \param[in] pmedata The PME structure
676 * \param[in] box The box matrix
677 * \param[in] stepWork Step schedule flags
678 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in the device memory.
679 * \param[in] lambdaQ The Coulomb lambda of the current state.
680 * \param[in] wcycle The wallcycle structure
682 static inline void launchPmeGpuSpread(gmx_pme_t* pmedata,
684 const StepWorkload& stepWork,
685 GpuEventSynchronizer* xReadyOnDevice,
687 gmx_wallcycle_t wcycle)
689 pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
690 pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle, lambdaQ);
693 /*! \brief Launch the FFT and gather stages of PME GPU
695 * This function only implements setting the output forces (no accumulation).
697 * \param[in] pmedata The PME structure
698 * \param[in] lambdaQ The Coulomb lambda of the current system state.
699 * \param[in] wcycle The wallcycle structure
700 * \param[in] stepWork Step schedule flags
702 static void launchPmeGpuFftAndGather(gmx_pme_t* pmedata,
704 gmx_wallcycle_t wcycle,
705 const gmx::StepWorkload& stepWork)
707 pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
708 pme_gpu_launch_gather(pmedata, wcycle, lambdaQ);
712 * Polling wait for either of the PME or nonbonded GPU tasks.
714 * Instead of a static order in waiting for GPU tasks, this function
715 * polls checking which of the two tasks completes first, and does the
716 * associated force buffer reduction overlapped with the other task.
717 * By doing that, unlike static scheduling order, it can always overlap
718 * one of the reductions, regardless of the GPU task completion order.
720 * \param[in] nbv Nonbonded verlet structure
721 * \param[in,out] pmedata PME module data
722 * \param[in,out] forceOutputs Output buffer for the forces and virial
723 * \param[in,out] enerd Energy data structure results are reduced into
724 * \param[in] lambdaQ The Coulomb lambda of the current system state.
725 * \param[in] stepWork Step schedule flags
726 * \param[in] wcycle The wallcycle structure
728 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
730 gmx::ForceOutputs* forceOutputs,
731 gmx_enerdata_t* enerd,
733 const StepWorkload& stepWork,
734 gmx_wallcycle_t wcycle)
736 bool isPmeGpuDone = false;
737 bool isNbGpuDone = false;
740 gmx::ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
741 gmx::ForceWithVirial& forceWithVirial = forceOutputs->forceWithVirial();
743 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
745 while (!isPmeGpuDone || !isNbGpuDone)
749 GpuTaskCompletion completionType =
750 (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
751 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, stepWork, wcycle, &forceWithVirial,
752 enerd, lambdaQ, completionType);
757 GpuTaskCompletion completionType =
758 (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
759 isNbGpuDone = Nbnxm::gpu_try_finish_task(
760 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
761 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(),
762 completionType, wcycle);
766 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShiftForces.force());
772 /*! \brief Set up the different force buffers; also does clearing.
774 * \param[in] forceHelperBuffers Helper force buffers
775 * \param[in] pull_work The pull work object.
776 * \param[in] inputrec input record
777 * \param[in] force force array
778 * \param[in] stepWork Step schedule flags
779 * \param[out] wcycle wallcycle recording structure
781 * \returns Cleared force output structure
783 static ForceOutputs setupForceOutputs(ForceHelperBuffers* forceHelperBuffers,
785 const t_inputrec& inputrec,
786 gmx::ArrayRefWithPadding<gmx::RVec> force,
787 const StepWorkload& stepWork,
788 gmx_wallcycle_t wcycle)
790 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
792 /* NOTE: We assume fr->shiftForces is all zeros here */
793 gmx::ForceWithShiftForces forceWithShiftForces(force, stepWork.computeVirial,
794 forceHelperBuffers->shiftForces());
796 if (stepWork.computeForces)
798 /* Clear the short- and long-range forces */
799 clearRVecs(forceWithShiftForces.force(), true);
801 /* Clear the shift forces */
802 clearRVecs(forceWithShiftForces.shiftForces(), false);
805 /* If we need to compute the virial, we might need a separate
806 * force buffer for algorithms for which the virial is calculated
807 * directly, such as PME. Otherwise, forceWithVirial uses the
808 * the same force (f in legacy calls) buffer as other algorithms.
810 const bool useSeparateForceWithVirialBuffer =
811 (stepWork.computeForces
812 && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
813 /* forceWithVirial uses the local atom range only */
814 gmx::ForceWithVirial forceWithVirial(
815 useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
816 : force.unpaddedArrayRef(),
817 stepWork.computeVirial);
819 if (useSeparateForceWithVirialBuffer)
821 /* TODO: update comment
822 * We only compute forces on local atoms. Note that vsites can
823 * spread to non-local atoms, but that part of the buffer is
824 * cleared separately in the vsite spreading code.
826 clearRVecs(forceWithVirial.force_, true);
829 if (inputrec.bPull && pull_have_constraint(pull_work))
831 clear_pull_forces(pull_work);
834 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
836 return ForceOutputs(forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(),
841 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
843 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec& inputrec,
844 const t_forcerec& fr,
845 const pull_t* pull_work,
847 const t_mdatoms& mdatoms,
848 const SimulationWorkload& simulationWork,
849 const StepWorkload& stepWork)
851 DomainLifetimeWorkload domainWork;
852 // Note that haveSpecialForces is constant over the whole run
853 domainWork.haveSpecialForces =
854 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
855 domainWork.haveCpuListedForceWork = false;
856 domainWork.haveCpuBondedWork = false;
857 for (const auto& listedForces : fr.listedForces)
859 if (listedForces.haveCpuListedForces(*fr.fcdata))
861 domainWork.haveCpuListedForceWork = true;
863 if (listedForces.haveCpuBondeds())
865 domainWork.haveCpuBondedWork = true;
868 domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
869 // Note that haveFreeEnergyWork is constant over the whole run
870 domainWork.haveFreeEnergyWork = (fr.efep != efepNO && mdatoms.nPerturbed != 0);
871 // We assume we have local force work if there are CPU
872 // force tasks including PME or nonbondeds.
873 domainWork.haveCpuLocalForceWork =
874 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
875 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
876 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
881 /*! \brief Set up force flag stuct from the force bitmask.
883 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
884 * \param[in] isNonbondedOn Global override, if false forces to turn off all nonbonded calculation.
885 * \param[in] simulationWork Simulation workload description.
886 * \param[in] rankHasPmeDuty If this rank computes PME.
888 * \returns New Stepworkload description.
890 static StepWorkload setupStepWorkload(const int legacyFlags,
891 const bool isNonbondedOn,
892 const SimulationWorkload& simulationWork,
893 const bool rankHasPmeDuty)
896 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
897 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
898 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
899 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
900 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
901 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
902 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
903 flags.computeNonbondedForces = ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && isNonbondedOn;
904 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
906 if (simulationWork.useGpuBufferOps)
908 GMX_ASSERT(simulationWork.useGpuNonbonded,
909 "Can only offload buffer ops if nonbonded computation is also offloaded");
911 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
912 // on virial steps the CPU reduction path is taken
913 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
914 flags.useGpuPmeFReduction = flags.useGpuFBufferOps
915 && (simulationWork.useGpuPme
916 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication));
922 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
924 * TODO: eliminate \p useGpuPmeOnThisRank when this is
925 * incorporated in DomainLifetimeWorkload.
927 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
928 gmx::GpuBonded* gpuBonded,
930 gmx_enerdata_t* enerd,
931 const gmx::MdrunScheduleWorkload& runScheduleWork,
932 bool useGpuPmeOnThisRank,
934 gmx_wallcycle_t wcycle)
936 if (runScheduleWork.simulationWork.useGpuNonbonded)
938 /* Launch pruning before buffer clearing because the API overhead of the
939 * clear kernel launches can leave the GPU idle while it could be running
942 if (nbv->isDynamicPruningStepGpu(step))
944 nbv->dispatchPruneKernelGpu(step);
947 /* now clear the GPU outputs while we finish the step on the CPU */
948 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
949 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
950 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
951 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
952 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
955 if (useGpuPmeOnThisRank)
957 pme_gpu_reinit_computation(pmedata, wcycle);
960 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
962 // in principle this should be included in the DD balancing region,
963 // but generally it is infrequent so we'll omit it for the sake of
965 gpuBonded->waitAccumulateEnergyTerms(enerd);
967 gpuBonded->clearEnergies();
971 //! \brief Data structure to hold dipole-related data and staging arrays
974 //! Dipole staging for fast summing over MPI
975 gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
976 //! Dipole staging for states A and B (index 0 and 1 resp.)
977 gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
981 static void reduceAndUpdateMuTot(DipoleData* dipoleData,
983 const bool haveFreeEnergy,
984 gmx::ArrayRef<const real> lambda,
986 const DDBalanceRegionHandler& ddBalanceRegionHandler)
990 gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
991 ddBalanceRegionHandler.reopenRegionCpu();
993 for (int i = 0; i < 2; i++)
995 for (int j = 0; j < DIM; j++)
997 dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
1001 if (!haveFreeEnergy)
1003 copy_rvec(dipoleData->muStateAB[0], muTotal);
1007 for (int j = 0; j < DIM; j++)
1009 muTotal[j] = (1.0 - lambda[efptCOUL]) * dipoleData->muStateAB[0][j]
1010 + lambda[efptCOUL] * dipoleData->muStateAB[1][j];
1015 void do_force(FILE* fplog,
1016 const t_commrec* cr,
1017 const gmx_multisim_t* ms,
1018 const t_inputrec* inputrec,
1020 gmx_enfrot* enforcedRotation,
1021 gmx::ImdSession* imdSession,
1025 gmx_wallcycle_t wcycle,
1026 const gmx_localtop_t* top,
1028 gmx::ArrayRefWithPadding<gmx::RVec> x,
1030 gmx::ForceBuffersView* forceView,
1032 const t_mdatoms* mdatoms,
1033 gmx_enerdata_t* enerd,
1034 gmx::ArrayRef<const real> lambda,
1036 gmx::MdrunScheduleWorkload* runScheduleWork,
1037 gmx::VirtualSitesHandler* vsite,
1042 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1044 auto force = forceView->forceWithPadding();
1045 GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1046 "The size of the force buffer should be at least the number of atoms to compute "
1049 nonbonded_verlet_t* nbv = fr->nbv.get();
1050 interaction_const_t* ic = fr->ic;
1051 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1053 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1056 runScheduleWork->stepWork = setupStepWorkload(legacyFlags, fr->bNonbonded, simulationWork,
1057 thisRankHasDuty(cr, DUTY_PME));
1058 const StepWorkload& stepWork = runScheduleWork->stepWork;
1061 const bool useGpuPmeOnThisRank = simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME);
1063 /* At a search step we need to start the first balancing region
1064 * somewhere early inside the step after communication during domain
1065 * decomposition (and not during the previous step as usual).
1067 if (stepWork.doNeighborSearch)
1069 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1072 clear_mat(vir_force);
1074 if (fr->pbcType != PbcType::No)
1076 /* Compute shift vectors every step,
1077 * because of pressure coupling or box deformation!
1079 if (stepWork.haveDynamicBox && stepWork.stateChanged)
1081 calc_shifts(box, fr->shift_vec);
1084 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1085 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1088 put_atoms_in_box_omp(fr->pbcType, box, x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1089 gmx_omp_nthreads_get(emntDefault));
1090 inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1094 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1096 const bool pmeSendCoordinatesFromGpu =
1097 GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1098 const bool reinitGpuPmePpComms =
1099 GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1101 const auto localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1102 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1103 AtomLocality::Local, simulationWork, stepWork)
1106 // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1107 // Otherwise the send will occur after H2D coordinate transfer.
1108 if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu)
1110 /* Send particle coordinates to the pme nodes */
1111 if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1113 GMX_RELEASE_ASSERT(false,
1114 "GPU update and separate PME ranks are only supported with GPU "
1115 "direct communication!");
1116 // TODO: when this code-path becomes supported add:
1117 // stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1120 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1121 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1122 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1123 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1126 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1127 // The local coordinates can be copied right away.
1128 // NOTE: Consider moving this copy to right after they are updated and constrained,
1129 // if the later is not offloaded.
1130 if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1132 if (stepWork.doNeighborSearch)
1134 // TODO refactor this to do_md, after partitioning.
1135 stateGpu->reinit(mdatoms->homenr,
1136 cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1137 if (useGpuPmeOnThisRank)
1139 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1140 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1143 // We need to copy coordinates when:
1144 // 1. Update is not offloaded
1145 // 2. The buffers were reinitialized on search step
1146 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1148 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1149 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1153 // TODO Update this comment when introducing SimulationWorkload
1155 // The conditions for gpuHaloExchange e.g. using GPU buffer
1156 // operations were checked before construction, so here we can
1157 // just use it and assert upon any conditions.
1158 const bool ddUsesGpuDirectCommunication =
1159 ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange.empty()));
1160 GMX_ASSERT(!ddUsesGpuDirectCommunication || stepWork.useGpuXBufferOps,
1161 "Must use coordinate buffer ops with GPU halo exchange");
1162 const bool useGpuForcesHaloExchange = ddUsesGpuDirectCommunication && stepWork.useGpuFBufferOps;
1164 // Copy coordinate from the GPU if update is on the GPU and there
1165 // are forces to be computed on the CPU, or for the computation of
1166 // virial, or if host-side data will be transferred from this task
1167 // to a remote task for halo exchange or PME-PP communication. At
1168 // search steps the current coordinates are already on the host,
1169 // hence copy is not needed.
1170 const bool haveHostPmePpComms =
1171 !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1172 const bool haveHostHaloExchangeComms = havePPDomainDecomposition(cr) && !ddUsesGpuDirectCommunication;
1174 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1175 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1176 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1177 || haveHostPmePpComms || haveHostHaloExchangeComms))
1179 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1180 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1181 haveCopiedXFromGpu = true;
1184 // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1185 // Otherwise the send will occur before the H2D coordinate transfer.
1186 if (!thisRankHasDuty(cr, DUTY_PME) && pmeSendCoordinatesFromGpu)
1188 /* Send particle coordinates to the pme nodes */
1189 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1190 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1191 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1192 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1195 if (useGpuPmeOnThisRank)
1197 launchPmeGpuSpread(fr->pmedata, box, stepWork, localXReadyOnDevice, lambda[efptCOUL], wcycle);
1200 /* do gridding for pair search */
1201 if (stepWork.doNeighborSearch)
1203 if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1205 fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1209 // - vzero is constant, do we need to pass it?
1210 // - box_diag should be passed directly to nbnxn_put_on_grid
1216 box_diag[XX] = box[XX][XX];
1217 box_diag[YY] = box[YY][YY];
1218 box_diag[ZZ] = box[ZZ][ZZ];
1220 wallcycle_start(wcycle, ewcNS);
1221 if (!DOMAINDECOMP(cr))
1223 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1224 nbnxn_put_on_grid(nbv, box, 0, vzero, box_diag, nullptr, { 0, mdatoms->homenr }, -1,
1225 fr->cginfo, x.unpaddedArrayRef(), 0, nullptr);
1226 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1230 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1231 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1232 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1235 nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1236 gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr), fr->cginfo);
1238 wallcycle_stop(wcycle, ewcNS);
1240 /* initialize the GPU nbnxm atom data and bonded data structures */
1241 if (simulationWork.useGpuNonbonded)
1243 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1245 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1246 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1247 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1251 /* Now we put all atoms on the grid, we can assign bonded
1252 * interactions to the GPU, where the grid order is
1253 * needed. Also the xq, f and fshift device buffers have
1254 * been reallocated if needed, so the bonded code can
1255 * learn about them. */
1256 // TODO the xq, f, and fshift buffers are now shared
1257 // resources, so they should be maintained by a
1258 // higher-level object than the nb module.
1259 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(
1260 nbv->getGridIndices(), top->idef, Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1261 Nbnxm::gpu_get_f(nbv->gpu_nbv), Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1263 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1266 // Need to run after the GPU-offload bonded interaction lists
1267 // are set up to be able to determine whether there is bonded work.
1268 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1269 *inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1271 wallcycle_start_nocount(wcycle, ewcNS);
1272 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1273 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1274 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1276 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1278 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1279 wallcycle_stop(wcycle, ewcNS);
1281 if (stepWork.useGpuXBufferOps)
1283 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1285 // For force buffer ops, we use the below conditon rather than
1286 // useGpuFBufferOps to ensure that init is performed even if this
1287 // NS step is also a virial step (on which f buf ops are deactivated).
1288 if (GMX_GPU_CUDA && simulationWork.useGpuBufferOps && simulationWork.useGpuNonbonded)
1290 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1291 nbv->atomdata_init_add_nbat_f_to_f_gpu(stateGpu->fReducedOnDevice());
1294 else if (!EI_TPI(inputrec->eI))
1296 if (stepWork.useGpuXBufferOps)
1298 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1299 nbv->convertCoordinatesGpu(AtomLocality::Local, false, stateGpu->getCoordinates(),
1300 localXReadyOnDevice);
1304 if (simulationWork.useGpuUpdate)
1306 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1307 GMX_ASSERT(haveCopiedXFromGpu,
1308 "a wait should only be triggered if copy has been scheduled");
1309 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1311 nbv->convertCoordinates(AtomLocality::Local, false, x.unpaddedArrayRef());
1315 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1317 if (simulationWork.useGpuNonbonded)
1319 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1321 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1323 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1324 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1325 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1327 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1329 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1330 // with X buffer ops offloaded to the GPU on all but the search steps
1332 // bonded work not split into separate local and non-local, so with DD
1333 // we can only launch the kernel after non-local coordinates have been received.
1334 if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1336 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1337 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1338 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1341 /* launch local nonbonded work on GPU */
1342 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1343 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1344 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1345 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1348 if (useGpuPmeOnThisRank)
1350 // In PME GPU and mixed mode we launch FFT / gather after the
1351 // X copy/transform to allow overlap as well as after the GPU NB
1352 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1353 // the nonbonded kernel.
1354 launchPmeGpuFftAndGather(fr->pmedata, lambda[efptCOUL], wcycle, stepWork);
1357 /* Communicate coordinates and sum dipole if necessary +
1358 do non-local pair search */
1359 if (havePPDomainDecomposition(cr))
1361 if (stepWork.doNeighborSearch)
1363 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1364 wallcycle_start_nocount(wcycle, ewcNS);
1365 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1366 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1367 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1369 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1370 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1371 wallcycle_stop(wcycle, ewcNS);
1372 // TODO refactor this GPU halo exchange re-initialisation
1373 // to location in do_md where GPU halo exchange is
1374 // constructed at partitioning, after above stateGpu
1375 // re-initialization has similarly been refactored
1376 if (ddUsesGpuDirectCommunication)
1378 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1383 if (ddUsesGpuDirectCommunication)
1385 // The following must be called after local setCoordinates (which records an event
1386 // when the coordinate data has been copied to the device).
1387 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1389 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1391 // non-local part of coordinate buffer must be copied back to host for CPU work
1392 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1397 // Note: GPU update + DD without direct communication is not supported,
1398 // a waitCoordinatesReadyOnHost() should be issued if it will be.
1399 GMX_ASSERT(!simulationWork.useGpuUpdate,
1400 "GPU update is not supported with CPU halo exchange");
1401 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1404 if (stepWork.useGpuXBufferOps)
1406 if (!useGpuPmeOnThisRank && !ddUsesGpuDirectCommunication)
1408 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1410 nbv->convertCoordinatesGpu(AtomLocality::NonLocal, false, stateGpu->getCoordinates(),
1411 stateGpu->getCoordinatesReadyOnDeviceEvent(
1412 AtomLocality::NonLocal, simulationWork, stepWork));
1416 nbv->convertCoordinates(AtomLocality::NonLocal, false, x.unpaddedArrayRef());
1420 if (simulationWork.useGpuNonbonded)
1422 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1424 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1426 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1427 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1428 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1431 if (domainWork.haveGpuBondedWork)
1433 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1434 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1435 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1438 /* launch non-local nonbonded tasks on GPU */
1439 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1440 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1442 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1444 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1448 if (simulationWork.useGpuNonbonded)
1450 /* launch D2H copy-back F */
1451 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1452 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1454 if (havePPDomainDecomposition(cr))
1456 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1458 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1459 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1461 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1463 fr->gpuBonded->launchEnergyTransfer();
1465 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1468 gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1469 if (fr->wholeMoleculeTransform)
1471 xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1474 DipoleData dipoleData;
1476 if (simulationWork.computeMuTot)
1478 const int start = 0;
1480 /* Calculate total (local) dipole moment in a temporary common array.
1481 * This makes it possible to sum them over nodes faster.
1483 gmx::ArrayRef<const gmx::RVec> xRef =
1484 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1485 calc_mu(start, mdatoms->homenr, xRef, mdatoms->chargeA, mdatoms->chargeB,
1486 mdatoms->nChargePerturbed, dipoleData.muStaging[0], dipoleData.muStaging[1]);
1488 reduceAndUpdateMuTot(&dipoleData, cr, (fr->efep != efepNO), lambda, muTotal, ddBalanceRegionHandler);
1491 /* Reset energies */
1492 reset_enerdata(enerd);
1494 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1496 wallcycle_start(wcycle, ewcPPDURINGPME);
1497 dd_force_flop_start(cr->dd, nrnb);
1500 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1501 // this wait ensures that the D2H transfer is complete.
1502 if ((simulationWork.useGpuUpdate)
1503 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1505 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1510 wallcycle_start(wcycle, ewcROT);
1511 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step,
1512 stepWork.doNeighborSearch);
1513 wallcycle_stop(wcycle, ewcROT);
1516 /* Start the force cycle counter.
1517 * Note that a different counter is used for dynamic load balancing.
1519 wallcycle_start(wcycle, ewcFORCE);
1521 // Set up and clear force outputs.
1522 // We use std::move to keep the compiler happy, it has no effect.
1523 ForceOutputs forceOut = setupForceOutputs(fr->forceHelperBuffers.get(), pull_work, *inputrec,
1524 std::move(force), stepWork, wcycle);
1526 /* We calculate the non-bonded forces, when done on the CPU, here.
1527 * We do this before calling do_force_lowlevel, because in that
1528 * function, the listed forces are calculated before PME, which
1529 * does communication. With this order, non-bonded and listed
1530 * force calculation imbalance can be balanced out by the domain
1531 * decomposition load balancing.
1534 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1536 if (!useOrEmulateGpuNb)
1538 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1541 if (fr->efep != efepNO)
1543 /* Calculate the local and non-local free energy interactions here.
1544 * Happens here on the CPU both with and without GPU.
1546 nbv->dispatchFreeEnergyKernel(InteractionLocality::Local, fr,
1547 as_rvec_array(x.unpaddedArrayRef().data()),
1548 &forceOut.forceWithShiftForces(), *mdatoms, inputrec->fepvals,
1549 lambda, enerd, stepWork, nrnb);
1551 if (havePPDomainDecomposition(cr))
1553 nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal, fr,
1554 as_rvec_array(x.unpaddedArrayRef().data()),
1555 &forceOut.forceWithShiftForces(), *mdatoms,
1556 inputrec->fepvals, lambda, enerd, stepWork, nrnb);
1560 if (!useOrEmulateGpuNb)
1562 if (havePPDomainDecomposition(cr))
1564 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1568 if (stepWork.computeForces)
1570 /* Add all the non-bonded force to the normal force array.
1571 * This can be split into a local and a non-local part when overlapping
1572 * communication with calculation with domain decomposition.
1574 wallcycle_stop(wcycle, ewcFORCE);
1575 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All, forceOut.forceWithShiftForces().force());
1576 wallcycle_start_nocount(wcycle, ewcFORCE);
1579 /* If there are multiple fshift output buffers we need to reduce them */
1580 if (stepWork.computeVirial)
1582 /* This is not in a subcounter because it takes a
1583 negligible and constant-sized amount of time */
1584 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat,
1585 forceOut.forceWithShiftForces().shiftForces());
1589 // TODO Force flags should include haveFreeEnergyWork for this domain
1590 if (ddUsesGpuDirectCommunication && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1592 /* Wait for non-local coordinate data to be copied from device */
1593 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1596 // Compute wall interactions, when present.
1597 // Note: should be moved to special forces.
1598 if (inputrec->nwall && stepWork.computeNonbondedForces)
1600 /* foreign lambda component for walls */
1601 real dvdl_walls = do_walls(*inputrec, *fr, box, *mdatoms, x.unpaddedConstArrayRef(),
1602 &forceOut.forceWithVirial(), lambda[efptVDW],
1603 enerd->grpp.ener[egLJSR].data(), nrnb);
1604 enerd->dvdl_lin[efptVDW] += dvdl_walls;
1607 if (stepWork.computeListedForces)
1609 /* Check whether we need to take into account PBC in listed interactions */
1610 bool needMolPbc = false;
1611 for (const auto& listedForces : fr->listedForces)
1613 if (listedForces.haveCpuListedForces(*fr->fcdata))
1615 needMolPbc = fr->bMolPBC;
1623 /* Since all atoms are in the rectangular or triclinic unit-cell,
1624 * only single box vector shifts (2 in x) are required.
1626 set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1629 for (auto& listedForces : fr->listedForces)
1631 listedForces.calculate(
1632 wcycle, box, inputrec->fepvals, cr, ms, x, xWholeMolecules, fr->fcdata.get(),
1633 hist, &forceOut, fr, &pbc, enerd, nrnb, lambda.data(), mdatoms,
1634 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr, stepWork);
1638 calculateLongRangeNonbondeds(fr, inputrec, cr, nrnb, wcycle, mdatoms, x.unpaddedConstArrayRef(),
1639 &forceOut.forceWithVirial(), enerd, box, lambda.data(),
1640 as_rvec_array(dipoleData.muStateAB), stepWork, ddBalanceRegionHandler);
1642 wallcycle_stop(wcycle, ewcFORCE);
1644 // VdW dispersion correction, only computed on master rank to avoid double counting
1645 if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1647 // Calculate long range corrections to pressure and energy
1648 const DispersionCorrection::Correction correction =
1649 fr->dispersionCorrection->calculate(box, lambda[efptVDW]);
1651 if (stepWork.computeEnergy)
1653 enerd->term[F_DISPCORR] = correction.energy;
1654 enerd->term[F_DVDL_VDW] += correction.dvdl;
1655 enerd->dvdl_lin[efptVDW] += correction.dvdl;
1657 if (stepWork.computeVirial)
1659 correction.correctVirial(vir_force);
1660 enerd->term[F_PDISPCORR] = correction.pressure;
1664 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation, imdSession, pull_work, step, t,
1665 wcycle, fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda,
1666 stepWork, &forceOut.forceWithVirial(), enerd, ed, stepWork.doNeighborSearch);
1669 // Will store the amount of cycles spent waiting for the GPU that
1670 // will be later used in the DLB accounting.
1671 float cycles_wait_gpu = 0;
1672 if (useOrEmulateGpuNb)
1674 auto& forceWithShiftForces = forceOut.forceWithShiftForces();
1676 /* wait for non-local forces (or calculate in emulation mode) */
1677 if (havePPDomainDecomposition(cr))
1679 if (simulationWork.useGpuNonbonded)
1681 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1682 nbv->gpu_nbv, stepWork, AtomLocality::NonLocal, enerd->grpp.ener[egLJSR].data(),
1683 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(), wcycle);
1687 wallcycle_start_nocount(wcycle, ewcFORCE);
1688 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes,
1689 step, nrnb, wcycle);
1690 wallcycle_stop(wcycle, ewcFORCE);
1693 if (stepWork.useGpuFBufferOps)
1695 gmx::FixedCapacityVector<GpuEventSynchronizer*, 1> dependencyList;
1697 // TODO: move this into DomainLifetimeWorkload, including the second part of the
1698 // condition The bonded and free energy CPU tasks can have non-local force
1699 // contributions which are a dependency for the GPU force reduction.
1700 bool haveNonLocalForceContribInCpuBuffer =
1701 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1703 if (haveNonLocalForceContribInCpuBuffer)
1705 stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(),
1706 AtomLocality::NonLocal);
1707 dependencyList.push_back(stateGpu->getForcesReadyOnDeviceEvent(
1708 AtomLocality::NonLocal, stepWork.useGpuFBufferOps));
1711 nbv->atomdata_add_nbat_f_to_f_gpu(AtomLocality::NonLocal, stateGpu->getForces(),
1712 pme_gpu_get_device_f(fr->pmedata), dependencyList,
1713 false, haveNonLocalForceContribInCpuBuffer);
1714 if (!useGpuForcesHaloExchange)
1716 // copy from GPU input for dd_move_f()
1717 stateGpu->copyForcesFromGpu(forceOut.forceWithShiftForces().force(),
1718 AtomLocality::NonLocal);
1723 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
1727 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
1729 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
1734 if (havePPDomainDecomposition(cr))
1736 /* We are done with the CPU compute.
1737 * We will now communicate the non-local forces.
1738 * If we use a GPU this will overlap with GPU work, so in that case
1739 * we do not close the DD force balancing region here.
1741 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1743 if (stepWork.computeForces)
1746 if (useGpuForcesHaloExchange)
1748 if (domainWork.haveCpuLocalForceWork)
1750 stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(), AtomLocality::Local);
1752 communicateGpuHaloForces(*cr, domainWork.haveCpuLocalForceWork);
1756 if (stepWork.useGpuFBufferOps)
1758 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
1760 dd_move_f(cr->dd, &forceOut.forceWithShiftForces(), wcycle);
1765 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1766 // an alternating wait/reduction scheme.
1767 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
1768 && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
1769 if (alternateGpuWait)
1771 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &forceOut, enerd, lambda[efptCOUL],
1775 if (!alternateGpuWait && useGpuPmeOnThisRank)
1777 pme_gpu_wait_and_reduce(fr->pmedata, stepWork, wcycle, &forceOut.forceWithVirial(), enerd,
1781 /* Wait for local GPU NB outputs on the non-alternating wait path */
1782 if (!alternateGpuWait && simulationWork.useGpuNonbonded)
1784 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1785 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1786 * but even with a step of 0.1 ms the difference is less than 1%
1789 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
1790 const float waitCycles = Nbnxm::gpu_wait_finish_task(
1791 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
1792 enerd->grpp.ener[egCOULSR].data(), forceOut.forceWithShiftForces().shiftForces(), wcycle);
1794 if (ddBalanceRegionHandler.useBalancingRegion())
1796 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1797 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
1799 /* We measured few cycles, it could be that the kernel
1800 * and transfer finished earlier and there was no actual
1801 * wait time, only API call overhead.
1802 * Then the actual time could be anywhere between 0 and
1803 * cycles_wait_est. We will use half of cycles_wait_est.
1805 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1807 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1811 if (fr->nbv->emulateGpu())
1813 // NOTE: emulation kernel is not included in the balancing region,
1814 // but emulation mode does not target performance anyway
1815 wallcycle_start_nocount(wcycle, ewcFORCE);
1816 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local,
1817 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes, step, nrnb, wcycle);
1818 wallcycle_stop(wcycle, ewcFORCE);
1821 // If on GPU PME-PP comms or GPU update path, receive forces from PME before GPU buffer ops
1822 // TODO refactor this and unify with below default-path call to the same function
1823 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME)
1824 && (simulationWork.useGpuPmePpCommunication || simulationWork.useGpuUpdate))
1826 /* In case of node-splitting, the PP nodes receive the long-range
1827 * forces, virial and energy from the PME nodes here.
1829 pme_receive_force_ener(fr, cr, &forceOut.forceWithVirial(), enerd,
1830 simulationWork.useGpuPmePpCommunication,
1831 stepWork.useGpuPmeFReduction, wcycle);
1835 /* Do the nonbonded GPU (or emulation) force buffer reduction
1836 * on the non-alternating path. */
1837 if (useOrEmulateGpuNb && !alternateGpuWait)
1839 // TODO simplify the below conditionals. Pass buffer and sync pointers at init stage rather than here. Unify getter fns for sameGPU/otherGPU cases.
1841 stepWork.useGpuPmeFReduction
1842 ? (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1843 : // PME force buffer on same GPU
1844 fr->pmePpCommGpu->getGpuForceStagingPtr()) // buffer received from other GPU
1845 : nullptr; // PME reduction not active on GPU
1847 GpuEventSynchronizer* const pmeSynchronizer =
1848 stepWork.useGpuPmeFReduction
1849 ? (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1850 : // PME force buffer on same GPU
1851 static_cast<GpuEventSynchronizer*>(
1852 fr->pmePpCommGpu->getForcesReadySynchronizer())) // buffer received from other GPU
1853 : nullptr; // PME reduction not active on GPU
1855 gmx::FixedCapacityVector<GpuEventSynchronizer*, 3> dependencyList;
1857 if (stepWork.useGpuPmeFReduction)
1859 dependencyList.push_back(pmeSynchronizer);
1862 gmx::ArrayRef<gmx::RVec> forceWithShift = forceOut.forceWithShiftForces().force();
1864 if (stepWork.useGpuFBufferOps)
1866 // Flag to specify whether the CPU force buffer has contributions to
1867 // local atoms. This depends on whether there are CPU-based force tasks
1868 // or when DD is active the halo exchange has resulted in contributions
1869 // from the non-local part.
1870 const bool haveLocalForceContribInCpuBuffer =
1871 (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
1873 // TODO: move these steps as early as possible:
1874 // - CPU f H2D should be as soon as all CPU-side forces are done
1875 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
1876 // before the next CPU task that consumes the forces: vsite spread or update)
1877 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
1878 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
1879 // These should be unified.
1880 if (haveLocalForceContribInCpuBuffer && !useGpuForcesHaloExchange)
1882 // Note: AtomLocality::All is used for the non-DD case because, as in this
1883 // case copyForcesToGpu() uses a separate stream, it allows overlap of
1884 // CPU force H2D with GPU force tasks on all streams including those in the
1885 // local stream which would otherwise be implicit dependencies for the
1886 // transfer and would not overlap.
1887 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
1889 stateGpu->copyForcesToGpu(forceWithShift, locality);
1890 dependencyList.push_back(
1891 stateGpu->getForcesReadyOnDeviceEvent(locality, stepWork.useGpuFBufferOps));
1893 if (useGpuForcesHaloExchange)
1895 dependencyList.push_back(cr->dd->gpuHaloExchange[0]->getForcesReadyOnDeviceEvent());
1897 nbv->atomdata_add_nbat_f_to_f_gpu(AtomLocality::Local, stateGpu->getForces(), pmeForcePtr,
1898 dependencyList, stepWork.useGpuPmeFReduction,
1899 haveLocalForceContribInCpuBuffer);
1900 // Copy forces to host if they are needed for update or if virtual sites are enabled.
1901 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
1902 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1903 // copy call done in sim_utils(...) for the output.
1904 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
1905 // they should not be copied in do_md(...) for the output.
1906 if (!simulationWork.useGpuUpdate || vsite)
1908 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
1909 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1914 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
1918 launchGpuEndOfStepTasks(nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork,
1919 useGpuPmeOnThisRank, step, wcycle);
1921 if (DOMAINDECOMP(cr))
1923 dd_force_flop_stop(cr->dd, nrnb);
1926 if (stepWork.computeForces)
1928 postProcessForceWithShiftForces(nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOut,
1929 vir_force, *mdatoms, *fr, vsite, stepWork);
1932 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
1933 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
1934 && !simulationWork.useGpuUpdate)
1936 /* In case of node-splitting, the PP nodes receive the long-range
1937 * forces, virial and energy from the PME nodes here.
1939 pme_receive_force_ener(fr, cr, &forceOut.forceWithVirial(), enerd,
1940 simulationWork.useGpuPmePpCommunication, false, wcycle);
1943 if (stepWork.computeForces)
1945 postProcessForces(cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOut, vir_force,
1946 mdatoms, fr, vsite, stepWork);
1949 if (stepWork.computeEnergy)
1951 /* Compute the final potential energy terms */
1952 accumulatePotentialEnergies(enerd, lambda, inputrec->fepvals);
1954 if (!EI_TPI(inputrec->eI))
1956 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1960 /* In case we don't have constraints and are using GPUs, the next balancing
1961 * region starts here.
1962 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1963 * virial calculation and COM pulling, is not thus not included in
1964 * the balance timing, which is ok as most tasks do communication.
1966 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);