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/applied_forces/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 */
407 /* When dynamic pair-list pruning is requested, we need to prune
408 * at nstlistPrune steps.
410 if (nbv->isDynamicPruningStepCpu(step))
412 /* Prune the pair-list beyond fr->ic->rlistPrune using
413 * the current coordinates of the atoms.
415 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
416 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
417 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
421 nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
424 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
426 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, v.ssize());
428 /* Note that we would like to avoid this conditional by putting it
429 * into the omp pragma instead, but then we still take the full
430 * omp parallel for overhead (at least with gcc5).
432 if (!useOpenmpThreading || nth == 1)
441 #pragma omp parallel for num_threads(nth) schedule(static)
442 for (gmx::index i = 0; i < v.ssize(); i++)
449 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
451 * \param groupOptions Group options, containing T-coupling options
453 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
455 real nrdfCoupled = 0;
456 real nrdfUncoupled = 0;
457 real kineticEnergy = 0;
458 for (int g = 0; g < groupOptions.ngtc; g++)
460 if (groupOptions.tau_t[g] >= 0)
462 nrdfCoupled += groupOptions.nrdf[g];
463 kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * BOLTZ;
467 nrdfUncoupled += groupOptions.nrdf[g];
471 /* This conditional with > also catches nrdf=0 */
472 if (nrdfCoupled > nrdfUncoupled)
474 return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
482 /*! \brief This routine checks that the potential energy is finite.
484 * Always checks that the potential energy is finite. If step equals
485 * inputrec.init_step also checks that the magnitude of the potential energy
486 * is reasonable. Terminates with a fatal error when a check fails.
487 * Note that passing this check does not guarantee finite forces,
488 * since those use slightly different arithmetics. But in most cases
489 * there is just a narrow coordinate range where forces are not finite
490 * and energies are finite.
492 * \param[in] step The step number, used for checking and printing
493 * \param[in] enerd The energy data; the non-bonded group energies need to be added to
494 * enerd.term[F_EPOT] before calling this routine \param[in] inputrec The input record
496 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
498 /* Threshold valid for comparing absolute potential energy against
499 * the kinetic energy. Normally one should not consider absolute
500 * potential energy values, but with a factor of one million
501 * we should never get false positives.
503 constexpr real c_thresholdFactor = 1e6;
505 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
506 real averageKineticEnergy = 0;
507 /* We only check for large potential energy at the initial step,
508 * because that is by far the most likely step for this too occur
509 * and because computing the average kinetic energy is not free.
510 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
511 * before they become NaN.
513 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
515 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
518 if (energyIsNotFinite
519 || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
524 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
525 "contributions to the energy are %g and %g, respectively. A %s potential energy "
526 "can be caused by overlapping interactions in bonded interactions or very large%s "
527 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
528 "configuration, incorrect interactions or parameters in the topology.",
529 step, enerd.term[F_EPOT], energyIsNotFinite ? "not finite" : "extremely high",
530 enerd.term[F_LJ], enerd.term[F_COUL_SR],
531 energyIsNotFinite ? "non-finite" : "very high", energyIsNotFinite ? " or Nan" : "");
535 /*! \brief Return true if there are special forces computed this step.
537 * The conditionals exactly correspond to those in computeSpecialForces().
539 static bool haveSpecialForces(const t_inputrec& inputrec,
540 const gmx::ForceProviders& forceProviders,
541 const pull_t* pull_work,
542 const bool computeForces,
546 return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
547 (inputrec.bPull && pull_have_potential(pull_work)) || // pull
548 inputrec.bRot || // enforced rotation
549 (ed != nullptr) || // flooding
550 (inputrec.bIMD && computeForces)); // IMD
553 /*! \brief Compute forces and/or energies for special algorithms
555 * The intention is to collect all calls to algorithms that compute
556 * forces on local atoms only and that do not contribute to the local
557 * virial sum (but add their virial contribution separately).
558 * Eventually these should likely all become ForceProviders.
559 * Within this function the intention is to have algorithms that do
560 * global communication at the end, so global barriers within the MD loop
561 * are as close together as possible.
563 * \param[in] fplog The log file
564 * \param[in] cr The communication record
565 * \param[in] inputrec The input record
566 * \param[in] awh The Awh module (nullptr if none in use).
567 * \param[in] enforcedRotation Enforced rotation module.
568 * \param[in] imdSession The IMD session
569 * \param[in] pull_work The pull work structure.
570 * \param[in] step The current MD step
571 * \param[in] t The current time
572 * \param[in,out] wcycle Wallcycle accounting struct
573 * \param[in,out] forceProviders Pointer to a list of force providers
574 * \param[in] box The unit cell
575 * \param[in] x The coordinates
576 * \param[in] mdatoms Per atom properties
577 * \param[in] lambda Array of free-energy lambda values
578 * \param[in] stepWork Step schedule flags
579 * \param[in,out] forceWithVirial Force and virial buffers
580 * \param[in,out] enerd Energy buffer
581 * \param[in,out] ed Essential dynamics pointer
582 * \param[in] didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
584 * \todo Remove didNeighborSearch, which is used incorrectly.
585 * \todo Convert all other algorithms called here to ForceProviders.
587 static void computeSpecialForces(FILE* fplog,
589 const t_inputrec* inputrec,
591 gmx_enfrot* enforcedRotation,
592 gmx::ImdSession* imdSession,
596 gmx_wallcycle_t wcycle,
597 gmx::ForceProviders* forceProviders,
599 gmx::ArrayRef<const gmx::RVec> x,
600 const t_mdatoms* mdatoms,
601 gmx::ArrayRef<const real> lambda,
602 const StepWorkload& stepWork,
603 gmx::ForceWithVirial* forceWithVirial,
604 gmx_enerdata_t* enerd,
606 bool didNeighborSearch)
608 /* NOTE: Currently all ForceProviders only provide forces.
609 * When they also provide energies, remove this conditional.
611 if (stepWork.computeForces)
613 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
614 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
616 /* Collect forces from modules */
617 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
620 if (inputrec->bPull && pull_have_potential(pull_work))
622 pull_potential_wrapper(cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work,
623 lambda.data(), t, wcycle);
627 const bool needForeignEnergyDifferences = awh->needForeignEnergyDifferences(step);
628 std::vector<double> foreignLambdaDeltaH, foreignLambdaDhDl;
629 if (needForeignEnergyDifferences)
631 enerd->foreignLambdaTerms.finalizePotentialContributions(enerd->dvdl_lin, lambda,
633 std::tie(foreignLambdaDeltaH, foreignLambdaDhDl) = enerd->foreignLambdaTerms.getTerms(cr);
636 enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(
637 inputrec->pbcType, mdatoms->massT, foreignLambdaDeltaH, foreignLambdaDhDl, box,
638 forceWithVirial, t, step, wcycle, fplog);
641 rvec* f = as_rvec_array(forceWithVirial->force_.data());
643 /* Add the forces from enforced rotation potentials (if any) */
646 wallcycle_start(wcycle, ewcROTadd);
647 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
648 wallcycle_stop(wcycle, ewcROTadd);
653 /* Note that since init_edsam() is called after the initialization
654 * of forcerec, edsam doesn't request the noVirSum force buffer.
655 * Thus if no other algorithm (e.g. PME) requires it, the forces
656 * here will contribute to the virial.
658 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
661 /* Add forces from interactive molecular dynamics (IMD), if any */
662 if (inputrec->bIMD && stepWork.computeForces)
664 imdSession->applyForces(f);
668 /*! \brief Launch the prepare_step and spread stages of PME GPU.
670 * \param[in] pmedata The PME structure
671 * \param[in] box The box matrix
672 * \param[in] stepWork Step schedule flags
673 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in the device memory.
674 * \param[in] lambdaQ The Coulomb lambda of the current state.
675 * \param[in] wcycle The wallcycle structure
677 static inline void launchPmeGpuSpread(gmx_pme_t* pmedata,
679 const StepWorkload& stepWork,
680 GpuEventSynchronizer* xReadyOnDevice,
682 gmx_wallcycle_t wcycle)
684 pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
685 pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle, lambdaQ);
688 /*! \brief Launch the FFT and gather stages of PME GPU
690 * This function only implements setting the output forces (no accumulation).
692 * \param[in] pmedata The PME structure
693 * \param[in] lambdaQ The Coulomb lambda of the current system state.
694 * \param[in] wcycle The wallcycle structure
695 * \param[in] stepWork Step schedule flags
697 static void launchPmeGpuFftAndGather(gmx_pme_t* pmedata,
699 gmx_wallcycle_t wcycle,
700 const gmx::StepWorkload& stepWork)
702 pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
703 pme_gpu_launch_gather(pmedata, wcycle, lambdaQ);
707 * Polling wait for either of the PME or nonbonded GPU tasks.
709 * Instead of a static order in waiting for GPU tasks, this function
710 * polls checking which of the two tasks completes first, and does the
711 * associated force buffer reduction overlapped with the other task.
712 * By doing that, unlike static scheduling order, it can always overlap
713 * one of the reductions, regardless of the GPU task completion order.
715 * \param[in] nbv Nonbonded verlet structure
716 * \param[in,out] pmedata PME module data
717 * \param[in,out] forceOutputs Output buffer for the forces and virial
718 * \param[in,out] enerd Energy data structure results are reduced into
719 * \param[in] lambdaQ The Coulomb lambda of the current system state.
720 * \param[in] stepWork Step schedule flags
721 * \param[in] wcycle The wallcycle structure
723 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
725 gmx::ForceOutputs* forceOutputs,
726 gmx_enerdata_t* enerd,
728 const StepWorkload& stepWork,
729 gmx_wallcycle_t wcycle)
731 bool isPmeGpuDone = false;
732 bool isNbGpuDone = false;
735 gmx::ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
736 gmx::ForceWithVirial& forceWithVirial = forceOutputs->forceWithVirial();
738 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
740 while (!isPmeGpuDone || !isNbGpuDone)
744 GpuTaskCompletion completionType =
745 (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
746 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, stepWork, wcycle, &forceWithVirial,
747 enerd, lambdaQ, completionType);
752 GpuTaskCompletion completionType =
753 (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
754 isNbGpuDone = Nbnxm::gpu_try_finish_task(
755 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
756 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(),
757 completionType, wcycle);
761 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShiftForces.force());
767 /*! \brief Set up the different force buffers; also does clearing.
769 * \param[in] forceHelperBuffers Helper force buffers
770 * \param[in] pull_work The pull work object.
771 * \param[in] inputrec input record
772 * \param[in] force force array
773 * \param[in] stepWork Step schedule flags
774 * \param[out] wcycle wallcycle recording structure
776 * \returns Cleared force output structure
778 static ForceOutputs setupForceOutputs(ForceHelperBuffers* forceHelperBuffers,
780 const t_inputrec& inputrec,
781 gmx::ArrayRefWithPadding<gmx::RVec> force,
782 const StepWorkload& stepWork,
783 gmx_wallcycle_t wcycle)
785 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
787 /* NOTE: We assume fr->shiftForces is all zeros here */
788 gmx::ForceWithShiftForces forceWithShiftForces(force, stepWork.computeVirial,
789 forceHelperBuffers->shiftForces());
791 if (stepWork.computeForces)
793 /* Clear the short- and long-range forces */
794 clearRVecs(forceWithShiftForces.force(), true);
796 /* Clear the shift forces */
797 clearRVecs(forceWithShiftForces.shiftForces(), false);
800 /* If we need to compute the virial, we might need a separate
801 * force buffer for algorithms for which the virial is calculated
802 * directly, such as PME. Otherwise, forceWithVirial uses the
803 * the same force (f in legacy calls) buffer as other algorithms.
805 const bool useSeparateForceWithVirialBuffer =
806 (stepWork.computeForces
807 && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
808 /* forceWithVirial uses the local atom range only */
809 gmx::ForceWithVirial forceWithVirial(
810 useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
811 : force.unpaddedArrayRef(),
812 stepWork.computeVirial);
814 if (useSeparateForceWithVirialBuffer)
816 /* TODO: update comment
817 * We only compute forces on local atoms. Note that vsites can
818 * spread to non-local atoms, but that part of the buffer is
819 * cleared separately in the vsite spreading code.
821 clearRVecs(forceWithVirial.force_, true);
824 if (inputrec.bPull && pull_have_constraint(pull_work))
826 clear_pull_forces(pull_work);
829 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
831 return ForceOutputs(forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(),
836 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
838 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec& inputrec,
839 const t_forcerec& fr,
840 const pull_t* pull_work,
842 const t_mdatoms& mdatoms,
843 const SimulationWorkload& simulationWork,
844 const StepWorkload& stepWork)
846 DomainLifetimeWorkload domainWork;
847 // Note that haveSpecialForces is constant over the whole run
848 domainWork.haveSpecialForces =
849 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
850 domainWork.haveCpuListedForceWork = false;
851 domainWork.haveCpuBondedWork = false;
852 for (const auto& listedForces : fr.listedForces)
854 if (listedForces.haveCpuListedForces(*fr.fcdata))
856 domainWork.haveCpuListedForceWork = true;
858 if (listedForces.haveCpuBondeds())
860 domainWork.haveCpuBondedWork = true;
863 domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
864 // Note that haveFreeEnergyWork is constant over the whole run
865 domainWork.haveFreeEnergyWork = (fr.efep != efepNO && mdatoms.nPerturbed != 0);
866 // We assume we have local force work if there are CPU
867 // force tasks including PME or nonbondeds.
868 domainWork.haveCpuLocalForceWork =
869 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
870 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
871 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
876 /*! \brief Set up force flag stuct from the force bitmask.
878 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
879 * \param[in] simulationWork Simulation workload description.
880 * \param[in] rankHasPmeDuty If this rank computes PME.
882 * \returns New Stepworkload description.
884 static StepWorkload setupStepWorkload(const int legacyFlags,
885 const SimulationWorkload& simulationWork,
886 const bool rankHasPmeDuty)
889 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
890 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
891 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
892 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
893 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
894 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
895 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
896 flags.computeNonbondedForces =
897 ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded;
898 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
900 if (simulationWork.useGpuBufferOps)
902 GMX_ASSERT(simulationWork.useGpuNonbonded,
903 "Can only offload buffer ops if nonbonded computation is also offloaded");
905 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
906 // on virial steps the CPU reduction path is taken
907 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
908 flags.useGpuPmeFReduction = flags.useGpuFBufferOps
909 && (simulationWork.useGpuPme
910 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication));
916 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
918 * TODO: eliminate \p useGpuPmeOnThisRank when this is
919 * incorporated in DomainLifetimeWorkload.
921 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
922 gmx::GpuBonded* gpuBonded,
924 gmx_enerdata_t* enerd,
925 const gmx::MdrunScheduleWorkload& runScheduleWork,
926 bool useGpuPmeOnThisRank,
928 gmx_wallcycle_t wcycle)
930 if (runScheduleWork.simulationWork.useGpuNonbonded)
932 /* Launch pruning before buffer clearing because the API overhead of the
933 * clear kernel launches can leave the GPU idle while it could be running
936 if (nbv->isDynamicPruningStepGpu(step))
938 nbv->dispatchPruneKernelGpu(step);
941 /* now clear the GPU outputs while we finish the step on the CPU */
942 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
943 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
944 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
945 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
946 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
949 if (useGpuPmeOnThisRank)
951 pme_gpu_reinit_computation(pmedata, wcycle);
954 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
956 // in principle this should be included in the DD balancing region,
957 // but generally it is infrequent so we'll omit it for the sake of
959 gpuBonded->waitAccumulateEnergyTerms(enerd);
961 gpuBonded->clearEnergies();
965 //! \brief Data structure to hold dipole-related data and staging arrays
968 //! Dipole staging for fast summing over MPI
969 gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
970 //! Dipole staging for states A and B (index 0 and 1 resp.)
971 gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
975 static void reduceAndUpdateMuTot(DipoleData* dipoleData,
977 const bool haveFreeEnergy,
978 gmx::ArrayRef<const real> lambda,
980 const DDBalanceRegionHandler& ddBalanceRegionHandler)
984 gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
985 ddBalanceRegionHandler.reopenRegionCpu();
987 for (int i = 0; i < 2; i++)
989 for (int j = 0; j < DIM; j++)
991 dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
997 copy_rvec(dipoleData->muStateAB[0], muTotal);
1001 for (int j = 0; j < DIM; j++)
1003 muTotal[j] = (1.0 - lambda[efptCOUL]) * dipoleData->muStateAB[0][j]
1004 + lambda[efptCOUL] * dipoleData->muStateAB[1][j];
1009 void do_force(FILE* fplog,
1010 const t_commrec* cr,
1011 const gmx_multisim_t* ms,
1012 const t_inputrec* inputrec,
1014 gmx_enfrot* enforcedRotation,
1015 gmx::ImdSession* imdSession,
1019 gmx_wallcycle_t wcycle,
1020 const gmx_localtop_t* top,
1022 gmx::ArrayRefWithPadding<gmx::RVec> x,
1024 gmx::ForceBuffersView* forceView,
1026 const t_mdatoms* mdatoms,
1027 gmx_enerdata_t* enerd,
1028 gmx::ArrayRef<const real> lambda,
1030 gmx::MdrunScheduleWorkload* runScheduleWork,
1031 gmx::VirtualSitesHandler* vsite,
1036 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1038 auto force = forceView->forceWithPadding();
1039 GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1040 "The size of the force buffer should be at least the number of atoms to compute "
1043 nonbonded_verlet_t* nbv = fr->nbv.get();
1044 interaction_const_t* ic = fr->ic;
1045 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1047 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1050 runScheduleWork->stepWork =
1051 setupStepWorkload(legacyFlags, simulationWork, thisRankHasDuty(cr, DUTY_PME));
1052 const StepWorkload& stepWork = runScheduleWork->stepWork;
1055 const bool useGpuPmeOnThisRank = simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME);
1057 /* At a search step we need to start the first balancing region
1058 * somewhere early inside the step after communication during domain
1059 * decomposition (and not during the previous step as usual).
1061 if (stepWork.doNeighborSearch)
1063 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1066 clear_mat(vir_force);
1068 if (fr->pbcType != PbcType::No)
1070 /* Compute shift vectors every step,
1071 * because of pressure coupling or box deformation!
1073 if (stepWork.haveDynamicBox && stepWork.stateChanged)
1075 calc_shifts(box, fr->shift_vec);
1078 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1079 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1082 put_atoms_in_box_omp(fr->pbcType, box, x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1083 gmx_omp_nthreads_get(emntDefault));
1084 inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1088 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1090 const bool pmeSendCoordinatesFromGpu =
1091 GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1092 const bool reinitGpuPmePpComms =
1093 GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1095 const auto localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1096 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1097 AtomLocality::Local, simulationWork, stepWork)
1100 // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1101 // Otherwise the send will occur after H2D coordinate transfer.
1102 if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu)
1104 /* Send particle coordinates to the pme nodes */
1105 if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1107 GMX_RELEASE_ASSERT(false,
1108 "GPU update and separate PME ranks are only supported with GPU "
1109 "direct communication!");
1110 // TODO: when this code-path becomes supported add:
1111 // stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1114 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1115 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1116 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1117 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1120 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1121 // The local coordinates can be copied right away.
1122 // NOTE: Consider moving this copy to right after they are updated and constrained,
1123 // if the later is not offloaded.
1124 if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1126 if (stepWork.doNeighborSearch)
1128 // TODO refactor this to do_md, after partitioning.
1129 stateGpu->reinit(mdatoms->homenr,
1130 cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1131 if (useGpuPmeOnThisRank)
1133 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1134 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1137 // We need to copy coordinates when:
1138 // 1. Update is not offloaded
1139 // 2. The buffers were reinitialized on search step
1140 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1142 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1143 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1147 // TODO Update this comment when introducing SimulationWorkload
1149 // The conditions for gpuHaloExchange e.g. using GPU buffer
1150 // operations were checked before construction, so here we can
1151 // just use it and assert upon any conditions.
1152 const bool ddUsesGpuDirectCommunication =
1153 ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange.empty()));
1154 GMX_ASSERT(!ddUsesGpuDirectCommunication || stepWork.useGpuXBufferOps,
1155 "Must use coordinate buffer ops with GPU halo exchange");
1156 const bool useGpuForcesHaloExchange = ddUsesGpuDirectCommunication && stepWork.useGpuFBufferOps;
1158 // Copy coordinate from the GPU if update is on the GPU and there
1159 // are forces to be computed on the CPU, or for the computation of
1160 // virial, or if host-side data will be transferred from this task
1161 // to a remote task for halo exchange or PME-PP communication. At
1162 // search steps the current coordinates are already on the host,
1163 // hence copy is not needed.
1164 const bool haveHostPmePpComms =
1165 !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1166 const bool haveHostHaloExchangeComms = havePPDomainDecomposition(cr) && !ddUsesGpuDirectCommunication;
1168 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1169 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1170 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1171 || haveHostPmePpComms || haveHostHaloExchangeComms))
1173 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1174 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1175 haveCopiedXFromGpu = true;
1178 // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1179 // Otherwise the send will occur before the H2D coordinate transfer.
1180 if (!thisRankHasDuty(cr, DUTY_PME) && pmeSendCoordinatesFromGpu)
1182 /* Send particle coordinates to the pme nodes */
1183 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1184 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1185 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1186 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1189 if (useGpuPmeOnThisRank)
1191 launchPmeGpuSpread(fr->pmedata, box, stepWork, localXReadyOnDevice, lambda[efptCOUL], wcycle);
1194 /* do gridding for pair search */
1195 if (stepWork.doNeighborSearch)
1197 if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1199 fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1203 // - vzero is constant, do we need to pass it?
1204 // - box_diag should be passed directly to nbnxn_put_on_grid
1210 box_diag[XX] = box[XX][XX];
1211 box_diag[YY] = box[YY][YY];
1212 box_diag[ZZ] = box[ZZ][ZZ];
1214 wallcycle_start(wcycle, ewcNS);
1215 if (!DOMAINDECOMP(cr))
1217 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1218 nbnxn_put_on_grid(nbv, box, 0, vzero, box_diag, nullptr, { 0, mdatoms->homenr }, -1,
1219 fr->cginfo, x.unpaddedArrayRef(), 0, nullptr);
1220 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1224 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1225 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1226 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1229 nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1230 gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr), fr->cginfo);
1232 wallcycle_stop(wcycle, ewcNS);
1234 /* initialize the GPU nbnxm atom data and bonded data structures */
1235 if (simulationWork.useGpuNonbonded)
1237 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1239 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1240 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1241 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1245 /* Now we put all atoms on the grid, we can assign bonded
1246 * interactions to the GPU, where the grid order is
1247 * needed. Also the xq, f and fshift device buffers have
1248 * been reallocated if needed, so the bonded code can
1249 * learn about them. */
1250 // TODO the xq, f, and fshift buffers are now shared
1251 // resources, so they should be maintained by a
1252 // higher-level object than the nb module.
1253 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(
1254 nbv->getGridIndices(), top->idef, Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1255 Nbnxm::gpu_get_f(nbv->gpu_nbv), Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1257 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1260 // Need to run after the GPU-offload bonded interaction lists
1261 // are set up to be able to determine whether there is bonded work.
1262 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1263 *inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1265 wallcycle_start_nocount(wcycle, ewcNS);
1266 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1267 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1268 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1270 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1272 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1273 wallcycle_stop(wcycle, ewcNS);
1275 if (stepWork.useGpuXBufferOps)
1277 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1279 // For force buffer ops, we use the below conditon rather than
1280 // useGpuFBufferOps to ensure that init is performed even if this
1281 // NS step is also a virial step (on which f buf ops are deactivated).
1282 if (GMX_GPU_CUDA && simulationWork.useGpuBufferOps && simulationWork.useGpuNonbonded)
1284 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1285 nbv->atomdata_init_add_nbat_f_to_f_gpu(stateGpu->fReducedOnDevice());
1288 else if (!EI_TPI(inputrec->eI))
1290 if (stepWork.useGpuXBufferOps)
1292 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1293 nbv->convertCoordinatesGpu(AtomLocality::Local, false, stateGpu->getCoordinates(),
1294 localXReadyOnDevice);
1298 if (simulationWork.useGpuUpdate)
1300 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1301 GMX_ASSERT(haveCopiedXFromGpu,
1302 "a wait should only be triggered if copy has been scheduled");
1303 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1305 nbv->convertCoordinates(AtomLocality::Local, false, x.unpaddedArrayRef());
1309 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1311 if (simulationWork.useGpuNonbonded)
1313 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1315 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1317 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1318 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1319 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1321 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1323 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1324 // with X buffer ops offloaded to the GPU on all but the search steps
1326 // bonded work not split into separate local and non-local, so with DD
1327 // we can only launch the kernel after non-local coordinates have been received.
1328 if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1330 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1331 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1332 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1335 /* launch local nonbonded work on GPU */
1336 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1337 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1338 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1339 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1342 if (useGpuPmeOnThisRank)
1344 // In PME GPU and mixed mode we launch FFT / gather after the
1345 // X copy/transform to allow overlap as well as after the GPU NB
1346 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1347 // the nonbonded kernel.
1348 launchPmeGpuFftAndGather(fr->pmedata, lambda[efptCOUL], wcycle, stepWork);
1351 /* Communicate coordinates and sum dipole if necessary +
1352 do non-local pair search */
1353 if (havePPDomainDecomposition(cr))
1355 if (stepWork.doNeighborSearch)
1357 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1358 wallcycle_start_nocount(wcycle, ewcNS);
1359 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1360 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1361 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1363 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1364 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1365 wallcycle_stop(wcycle, ewcNS);
1366 // TODO refactor this GPU halo exchange re-initialisation
1367 // to location in do_md where GPU halo exchange is
1368 // constructed at partitioning, after above stateGpu
1369 // re-initialization has similarly been refactored
1370 if (ddUsesGpuDirectCommunication)
1372 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1377 if (ddUsesGpuDirectCommunication)
1379 // The following must be called after local setCoordinates (which records an event
1380 // when the coordinate data has been copied to the device).
1381 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1383 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1385 // non-local part of coordinate buffer must be copied back to host for CPU work
1386 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1391 // Note: GPU update + DD without direct communication is not supported,
1392 // a waitCoordinatesReadyOnHost() should be issued if it will be.
1393 GMX_ASSERT(!simulationWork.useGpuUpdate,
1394 "GPU update is not supported with CPU halo exchange");
1395 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1398 if (stepWork.useGpuXBufferOps)
1400 if (!useGpuPmeOnThisRank && !ddUsesGpuDirectCommunication)
1402 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1404 nbv->convertCoordinatesGpu(AtomLocality::NonLocal, false, stateGpu->getCoordinates(),
1405 stateGpu->getCoordinatesReadyOnDeviceEvent(
1406 AtomLocality::NonLocal, simulationWork, stepWork));
1410 nbv->convertCoordinates(AtomLocality::NonLocal, false, x.unpaddedArrayRef());
1414 if (simulationWork.useGpuNonbonded)
1416 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1418 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1420 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1421 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1422 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1425 if (domainWork.haveGpuBondedWork)
1427 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1428 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1429 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1432 /* launch non-local nonbonded tasks on GPU */
1433 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1434 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1436 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1438 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1442 if (simulationWork.useGpuNonbonded)
1444 /* launch D2H copy-back F */
1445 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1446 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1448 if (havePPDomainDecomposition(cr))
1450 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1452 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1453 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1455 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1457 fr->gpuBonded->launchEnergyTransfer();
1459 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1462 gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1463 if (fr->wholeMoleculeTransform)
1465 xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1468 DipoleData dipoleData;
1470 if (simulationWork.computeMuTot)
1472 const int start = 0;
1474 /* Calculate total (local) dipole moment in a temporary common array.
1475 * This makes it possible to sum them over nodes faster.
1477 gmx::ArrayRef<const gmx::RVec> xRef =
1478 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1479 calc_mu(start, mdatoms->homenr, xRef, mdatoms->chargeA, mdatoms->chargeB,
1480 mdatoms->nChargePerturbed, dipoleData.muStaging[0], dipoleData.muStaging[1]);
1482 reduceAndUpdateMuTot(&dipoleData, cr, (fr->efep != efepNO), lambda, muTotal, ddBalanceRegionHandler);
1485 /* Reset energies */
1486 reset_enerdata(enerd);
1488 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1490 wallcycle_start(wcycle, ewcPPDURINGPME);
1491 dd_force_flop_start(cr->dd, nrnb);
1494 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1495 // this wait ensures that the D2H transfer is complete.
1496 if ((simulationWork.useGpuUpdate)
1497 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1499 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1504 wallcycle_start(wcycle, ewcROT);
1505 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step,
1506 stepWork.doNeighborSearch);
1507 wallcycle_stop(wcycle, ewcROT);
1510 /* Start the force cycle counter.
1511 * Note that a different counter is used for dynamic load balancing.
1513 wallcycle_start(wcycle, ewcFORCE);
1515 // Set up and clear force outputs.
1516 // We use std::move to keep the compiler happy, it has no effect.
1517 ForceOutputs forceOut = setupForceOutputs(fr->forceHelperBuffers.get(), pull_work, *inputrec,
1518 std::move(force), stepWork, wcycle);
1520 /* We calculate the non-bonded forces, when done on the CPU, here.
1521 * We do this before calling do_force_lowlevel, because in that
1522 * function, the listed forces are calculated before PME, which
1523 * does communication. With this order, non-bonded and listed
1524 * force calculation imbalance can be balanced out by the domain
1525 * decomposition load balancing.
1528 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1530 if (!useOrEmulateGpuNb)
1532 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1535 if (fr->efep != efepNO)
1537 /* Calculate the local and non-local free energy interactions here.
1538 * Happens here on the CPU both with and without GPU.
1540 nbv->dispatchFreeEnergyKernel(InteractionLocality::Local, fr,
1541 as_rvec_array(x.unpaddedArrayRef().data()),
1542 &forceOut.forceWithShiftForces(), *mdatoms, inputrec->fepvals,
1543 lambda, enerd, stepWork, nrnb);
1545 if (havePPDomainDecomposition(cr))
1547 nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal, fr,
1548 as_rvec_array(x.unpaddedArrayRef().data()),
1549 &forceOut.forceWithShiftForces(), *mdatoms,
1550 inputrec->fepvals, lambda, enerd, stepWork, nrnb);
1554 if (!useOrEmulateGpuNb)
1556 if (havePPDomainDecomposition(cr))
1558 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1562 if (stepWork.computeForces)
1564 /* Add all the non-bonded force to the normal force array.
1565 * This can be split into a local and a non-local part when overlapping
1566 * communication with calculation with domain decomposition.
1568 wallcycle_stop(wcycle, ewcFORCE);
1569 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All, forceOut.forceWithShiftForces().force());
1570 wallcycle_start_nocount(wcycle, ewcFORCE);
1573 /* If there are multiple fshift output buffers we need to reduce them */
1574 if (stepWork.computeVirial)
1576 /* This is not in a subcounter because it takes a
1577 negligible and constant-sized amount of time */
1578 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat,
1579 forceOut.forceWithShiftForces().shiftForces());
1583 // TODO Force flags should include haveFreeEnergyWork for this domain
1584 if (ddUsesGpuDirectCommunication && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1586 /* Wait for non-local coordinate data to be copied from device */
1587 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1590 // Compute wall interactions, when present.
1591 // Note: should be moved to special forces.
1592 if (inputrec->nwall && stepWork.computeNonbondedForces)
1594 /* foreign lambda component for walls */
1595 real dvdl_walls = do_walls(*inputrec, *fr, box, *mdatoms, x.unpaddedConstArrayRef(),
1596 &forceOut.forceWithVirial(), lambda[efptVDW],
1597 enerd->grpp.ener[egLJSR].data(), nrnb);
1598 enerd->dvdl_lin[efptVDW] += dvdl_walls;
1601 if (stepWork.computeListedForces)
1603 /* Check whether we need to take into account PBC in listed interactions */
1604 bool needMolPbc = false;
1605 for (const auto& listedForces : fr->listedForces)
1607 if (listedForces.haveCpuListedForces(*fr->fcdata))
1609 needMolPbc = fr->bMolPBC;
1617 /* Since all atoms are in the rectangular or triclinic unit-cell,
1618 * only single box vector shifts (2 in x) are required.
1620 set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1623 for (auto& listedForces : fr->listedForces)
1625 listedForces.calculate(
1626 wcycle, box, inputrec->fepvals, cr, ms, x, xWholeMolecules, fr->fcdata.get(),
1627 hist, &forceOut, fr, &pbc, enerd, nrnb, lambda.data(), mdatoms,
1628 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr, stepWork);
1632 calculateLongRangeNonbondeds(fr, inputrec, cr, nrnb, wcycle, mdatoms, x.unpaddedConstArrayRef(),
1633 &forceOut.forceWithVirial(), enerd, box, lambda.data(),
1634 as_rvec_array(dipoleData.muStateAB), stepWork, ddBalanceRegionHandler);
1636 wallcycle_stop(wcycle, ewcFORCE);
1638 // VdW dispersion correction, only computed on master rank to avoid double counting
1639 if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1641 // Calculate long range corrections to pressure and energy
1642 const DispersionCorrection::Correction correction =
1643 fr->dispersionCorrection->calculate(box, lambda[efptVDW]);
1645 if (stepWork.computeEnergy)
1647 enerd->term[F_DISPCORR] = correction.energy;
1648 enerd->term[F_DVDL_VDW] += correction.dvdl;
1649 enerd->dvdl_lin[efptVDW] += correction.dvdl;
1651 if (stepWork.computeVirial)
1653 correction.correctVirial(vir_force);
1654 enerd->term[F_PDISPCORR] = correction.pressure;
1658 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation, imdSession, pull_work, step, t,
1659 wcycle, fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda,
1660 stepWork, &forceOut.forceWithVirial(), enerd, ed, stepWork.doNeighborSearch);
1663 // Will store the amount of cycles spent waiting for the GPU that
1664 // will be later used in the DLB accounting.
1665 float cycles_wait_gpu = 0;
1666 if (useOrEmulateGpuNb)
1668 auto& forceWithShiftForces = forceOut.forceWithShiftForces();
1670 /* wait for non-local forces (or calculate in emulation mode) */
1671 if (havePPDomainDecomposition(cr))
1673 if (simulationWork.useGpuNonbonded)
1675 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1676 nbv->gpu_nbv, stepWork, AtomLocality::NonLocal, enerd->grpp.ener[egLJSR].data(),
1677 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(), wcycle);
1681 wallcycle_start_nocount(wcycle, ewcFORCE);
1682 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes,
1683 step, nrnb, wcycle);
1684 wallcycle_stop(wcycle, ewcFORCE);
1687 if (stepWork.useGpuFBufferOps)
1689 gmx::FixedCapacityVector<GpuEventSynchronizer*, 1> dependencyList;
1691 // TODO: move this into DomainLifetimeWorkload, including the second part of the
1692 // condition The bonded and free energy CPU tasks can have non-local force
1693 // contributions which are a dependency for the GPU force reduction.
1694 bool haveNonLocalForceContribInCpuBuffer =
1695 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1697 if (haveNonLocalForceContribInCpuBuffer)
1699 stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(),
1700 AtomLocality::NonLocal);
1701 dependencyList.push_back(stateGpu->getForcesReadyOnDeviceEvent(
1702 AtomLocality::NonLocal, stepWork.useGpuFBufferOps));
1705 nbv->atomdata_add_nbat_f_to_f_gpu(AtomLocality::NonLocal, stateGpu->getForces(),
1706 pme_gpu_get_device_f(fr->pmedata), dependencyList,
1707 false, haveNonLocalForceContribInCpuBuffer);
1708 if (!useGpuForcesHaloExchange)
1710 // copy from GPU input for dd_move_f()
1711 stateGpu->copyForcesFromGpu(forceOut.forceWithShiftForces().force(),
1712 AtomLocality::NonLocal);
1717 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
1721 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
1723 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
1728 if (havePPDomainDecomposition(cr))
1730 /* We are done with the CPU compute.
1731 * We will now communicate the non-local forces.
1732 * If we use a GPU this will overlap with GPU work, so in that case
1733 * we do not close the DD force balancing region here.
1735 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1737 if (stepWork.computeForces)
1740 if (useGpuForcesHaloExchange)
1742 if (domainWork.haveCpuLocalForceWork)
1744 stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(), AtomLocality::Local);
1746 communicateGpuHaloForces(*cr, domainWork.haveCpuLocalForceWork);
1750 if (stepWork.useGpuFBufferOps)
1752 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
1754 dd_move_f(cr->dd, &forceOut.forceWithShiftForces(), wcycle);
1759 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1760 // an alternating wait/reduction scheme.
1761 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
1762 && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
1763 if (alternateGpuWait)
1765 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &forceOut, enerd, lambda[efptCOUL],
1769 if (!alternateGpuWait && useGpuPmeOnThisRank)
1771 pme_gpu_wait_and_reduce(fr->pmedata, stepWork, wcycle, &forceOut.forceWithVirial(), enerd,
1775 /* Wait for local GPU NB outputs on the non-alternating wait path */
1776 if (!alternateGpuWait && simulationWork.useGpuNonbonded)
1778 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1779 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1780 * but even with a step of 0.1 ms the difference is less than 1%
1783 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
1784 const float waitCycles = Nbnxm::gpu_wait_finish_task(
1785 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
1786 enerd->grpp.ener[egCOULSR].data(), forceOut.forceWithShiftForces().shiftForces(), wcycle);
1788 if (ddBalanceRegionHandler.useBalancingRegion())
1790 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1791 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
1793 /* We measured few cycles, it could be that the kernel
1794 * and transfer finished earlier and there was no actual
1795 * wait time, only API call overhead.
1796 * Then the actual time could be anywhere between 0 and
1797 * cycles_wait_est. We will use half of cycles_wait_est.
1799 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1801 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1805 if (fr->nbv->emulateGpu())
1807 // NOTE: emulation kernel is not included in the balancing region,
1808 // but emulation mode does not target performance anyway
1809 wallcycle_start_nocount(wcycle, ewcFORCE);
1810 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local,
1811 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes, step, nrnb, wcycle);
1812 wallcycle_stop(wcycle, ewcFORCE);
1815 // If on GPU PME-PP comms or GPU update path, receive forces from PME before GPU buffer ops
1816 // TODO refactor this and unify with below default-path call to the same function
1817 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME)
1818 && (simulationWork.useGpuPmePpCommunication || simulationWork.useGpuUpdate))
1820 /* In case of node-splitting, the PP nodes receive the long-range
1821 * forces, virial and energy from the PME nodes here.
1823 pme_receive_force_ener(fr, cr, &forceOut.forceWithVirial(), enerd,
1824 simulationWork.useGpuPmePpCommunication,
1825 stepWork.useGpuPmeFReduction, wcycle);
1829 /* Do the nonbonded GPU (or emulation) force buffer reduction
1830 * on the non-alternating path. */
1831 if (useOrEmulateGpuNb && !alternateGpuWait)
1833 // TODO simplify the below conditionals. Pass buffer and sync pointers at init stage rather than here. Unify getter fns for sameGPU/otherGPU cases.
1835 stepWork.useGpuPmeFReduction
1836 ? (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1837 : // PME force buffer on same GPU
1838 fr->pmePpCommGpu->getGpuForceStagingPtr()) // buffer received from other GPU
1839 : nullptr; // PME reduction not active on GPU
1841 GpuEventSynchronizer* const pmeSynchronizer =
1842 stepWork.useGpuPmeFReduction
1843 ? (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1844 : // PME force buffer on same GPU
1845 static_cast<GpuEventSynchronizer*>(
1846 fr->pmePpCommGpu->getForcesReadySynchronizer())) // buffer received from other GPU
1847 : nullptr; // PME reduction not active on GPU
1849 gmx::FixedCapacityVector<GpuEventSynchronizer*, 3> dependencyList;
1851 if (stepWork.useGpuPmeFReduction)
1853 dependencyList.push_back(pmeSynchronizer);
1856 gmx::ArrayRef<gmx::RVec> forceWithShift = forceOut.forceWithShiftForces().force();
1858 if (stepWork.useGpuFBufferOps)
1860 // Flag to specify whether the CPU force buffer has contributions to
1861 // local atoms. This depends on whether there are CPU-based force tasks
1862 // or when DD is active the halo exchange has resulted in contributions
1863 // from the non-local part.
1864 const bool haveLocalForceContribInCpuBuffer =
1865 (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
1867 // TODO: move these steps as early as possible:
1868 // - CPU f H2D should be as soon as all CPU-side forces are done
1869 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
1870 // before the next CPU task that consumes the forces: vsite spread or update)
1871 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
1872 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
1873 // These should be unified.
1874 if (haveLocalForceContribInCpuBuffer && !useGpuForcesHaloExchange)
1876 // Note: AtomLocality::All is used for the non-DD case because, as in this
1877 // case copyForcesToGpu() uses a separate stream, it allows overlap of
1878 // CPU force H2D with GPU force tasks on all streams including those in the
1879 // local stream which would otherwise be implicit dependencies for the
1880 // transfer and would not overlap.
1881 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
1883 stateGpu->copyForcesToGpu(forceWithShift, locality);
1884 dependencyList.push_back(
1885 stateGpu->getForcesReadyOnDeviceEvent(locality, stepWork.useGpuFBufferOps));
1887 if (useGpuForcesHaloExchange)
1889 dependencyList.push_back(cr->dd->gpuHaloExchange[0]->getForcesReadyOnDeviceEvent());
1891 nbv->atomdata_add_nbat_f_to_f_gpu(AtomLocality::Local, stateGpu->getForces(), pmeForcePtr,
1892 dependencyList, stepWork.useGpuPmeFReduction,
1893 haveLocalForceContribInCpuBuffer);
1894 // Copy forces to host if they are needed for update or if virtual sites are enabled.
1895 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
1896 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1897 // copy call done in sim_utils(...) for the output.
1898 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
1899 // they should not be copied in do_md(...) for the output.
1900 if (!simulationWork.useGpuUpdate || vsite)
1902 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
1903 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1908 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
1912 launchGpuEndOfStepTasks(nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork,
1913 useGpuPmeOnThisRank, step, wcycle);
1915 if (DOMAINDECOMP(cr))
1917 dd_force_flop_stop(cr->dd, nrnb);
1920 if (stepWork.computeForces)
1922 postProcessForceWithShiftForces(nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOut,
1923 vir_force, *mdatoms, *fr, vsite, stepWork);
1926 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
1927 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
1928 && !simulationWork.useGpuUpdate)
1930 /* In case of node-splitting, the PP nodes receive the long-range
1931 * forces, virial and energy from the PME nodes here.
1933 pme_receive_force_ener(fr, cr, &forceOut.forceWithVirial(), enerd,
1934 simulationWork.useGpuPmePpCommunication, false, wcycle);
1937 if (stepWork.computeForces)
1939 postProcessForces(cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOut, vir_force,
1940 mdatoms, fr, vsite, stepWork);
1943 if (stepWork.computeEnergy)
1945 /* Compute the final potential energy terms */
1946 accumulatePotentialEnergies(enerd, lambda, inputrec->fepvals);
1948 if (!EI_TPI(inputrec->eI))
1950 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1954 /* In case we don't have constraints and are using GPUs, the next balancing
1955 * region starts here.
1956 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1957 * virial calculation and COM pulling, is not thus not included in
1958 * the balance timing, which is ok as most tasks do communication.
1960 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);