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.
49 #include "gromacs/applied_forces/awh/awh.h"
50 #include "gromacs/domdec/dlbtiming.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/domdec/gpuhaloexchange.h"
54 #include "gromacs/domdec/partition.h"
55 #include "gromacs/essentialdynamics/edsam.h"
56 #include "gromacs/ewald/pme.h"
57 #include "gromacs/ewald/pme_pp.h"
58 #include "gromacs/ewald/pme_pp_comm_gpu.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
61 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
62 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
63 #include "gromacs/gpu_utils/gpu_utils.h"
64 #include "gromacs/imd/imd.h"
65 #include "gromacs/listed_forces/disre.h"
66 #include "gromacs/listed_forces/gpubonded.h"
67 #include "gromacs/listed_forces/listed_forces.h"
68 #include "gromacs/listed_forces/orires.h"
69 #include "gromacs/math/arrayrefwithpadding.h"
70 #include "gromacs/math/functions.h"
71 #include "gromacs/math/units.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/math/vecdump.h"
74 #include "gromacs/mdlib/calcmu.h"
75 #include "gromacs/mdlib/calcvir.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/dispersioncorrection.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/force.h"
80 #include "gromacs/mdlib/force_flags.h"
81 #include "gromacs/mdlib/forcerec.h"
82 #include "gromacs/mdlib/gmx_omp_nthreads.h"
83 #include "gromacs/mdlib/update.h"
84 #include "gromacs/mdlib/vsite.h"
85 #include "gromacs/mdlib/wall.h"
86 #include "gromacs/mdlib/wholemoleculetransform.h"
87 #include "gromacs/mdtypes/commrec.h"
88 #include "gromacs/mdtypes/enerdata.h"
89 #include "gromacs/mdtypes/forcebuffers.h"
90 #include "gromacs/mdtypes/forceoutput.h"
91 #include "gromacs/mdtypes/forcerec.h"
92 #include "gromacs/mdtypes/iforceprovider.h"
93 #include "gromacs/mdtypes/inputrec.h"
94 #include "gromacs/mdtypes/md_enums.h"
95 #include "gromacs/mdtypes/mdatom.h"
96 #include "gromacs/mdtypes/multipletimestepping.h"
97 #include "gromacs/mdtypes/simulation_workload.h"
98 #include "gromacs/mdtypes/state.h"
99 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
100 #include "gromacs/nbnxm/gpu_data_mgmt.h"
101 #include "gromacs/nbnxm/nbnxm.h"
102 #include "gromacs/nbnxm/nbnxm_gpu.h"
103 #include "gromacs/pbcutil/ishift.h"
104 #include "gromacs/pbcutil/pbc.h"
105 #include "gromacs/pulling/pull.h"
106 #include "gromacs/pulling/pull_rotation.h"
107 #include "gromacs/timing/cyclecounter.h"
108 #include "gromacs/timing/gpu_timing.h"
109 #include "gromacs/timing/wallcycle.h"
110 #include "gromacs/timing/wallcyclereporting.h"
111 #include "gromacs/timing/walltime_accounting.h"
112 #include "gromacs/topology/topology.h"
113 #include "gromacs/utility/arrayref.h"
114 #include "gromacs/utility/basedefinitions.h"
115 #include "gromacs/utility/cstringutil.h"
116 #include "gromacs/utility/exceptions.h"
117 #include "gromacs/utility/fatalerror.h"
118 #include "gromacs/utility/fixedcapacityvector.h"
119 #include "gromacs/utility/gmxassert.h"
120 #include "gromacs/utility/gmxmpi.h"
121 #include "gromacs/utility/logger.h"
122 #include "gromacs/utility/smalloc.h"
123 #include "gromacs/utility/strconvert.h"
124 #include "gromacs/utility/sysinfo.h"
126 #include "gpuforcereduction.h"
129 using gmx::AtomLocality;
130 using gmx::DomainLifetimeWorkload;
131 using gmx::ForceOutputs;
132 using gmx::ForceWithShiftForces;
133 using gmx::InteractionLocality;
135 using gmx::SimulationWorkload;
136 using gmx::StepWorkload;
138 // TODO: this environment variable allows us to verify before release
139 // that on less common architectures the total cost of polling is not larger than
140 // a blocking wait (so polling does not introduce overhead when the static
141 // PME-first ordering would suffice).
142 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
144 static void sum_forces(ArrayRef<RVec> f, ArrayRef<const RVec> forceToAdd)
146 GMX_ASSERT(f.size() >= forceToAdd.size(), "Accumulation buffer should be sufficiently large");
147 const int end = forceToAdd.size();
149 int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
150 #pragma omp parallel for num_threads(nt) schedule(static)
151 for (int i = 0; i < end; i++)
153 rvec_inc(f[i], forceToAdd[i]);
157 static void calc_virial(int start,
160 const gmx::ForceWithShiftForces& forceWithShiftForces,
164 const t_forcerec* fr,
167 /* The short-range virial from surrounding boxes */
168 const rvec* fshift = as_rvec_array(forceWithShiftForces.shiftForces().data());
169 calc_vir(SHIFTS, fr->shift_vec, fshift, vir_part, pbcType == PbcType::Screw, box);
170 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
172 /* Calculate partial virial, for local atoms only, based on short range.
173 * Total virial is computed in global_stat, called from do_md
175 const rvec* f = as_rvec_array(forceWithShiftForces.force().data());
176 f_calc_vir(start, start + homenr, x, f, vir_part, box);
177 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
181 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
185 static void pull_potential_wrapper(const t_commrec* cr,
186 const t_inputrec* ir,
188 gmx::ArrayRef<const gmx::RVec> x,
189 gmx::ForceWithVirial* force,
190 const t_mdatoms* mdatoms,
191 gmx_enerdata_t* enerd,
195 gmx_wallcycle_t wcycle)
200 /* Calculate the center of mass forces, this requires communication,
201 * which is why pull_potential is called close to other communication.
203 wallcycle_start(wcycle, ewcPULLPOT);
204 set_pbc(&pbc, ir->pbcType, box);
206 enerd->term[F_COM_PULL] +=
207 pull_potential(pull_work, mdatoms->massT, &pbc, cr, t, lambda[efptRESTRAINT],
208 as_rvec_array(x.data()), force, &dvdl);
209 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
210 wallcycle_stop(wcycle, ewcPULLPOT);
213 static void pme_receive_force_ener(t_forcerec* fr,
215 gmx::ForceWithVirial* forceWithVirial,
216 gmx_enerdata_t* enerd,
217 bool useGpuPmePpComms,
218 bool receivePmeForceToGpu,
219 gmx_wallcycle_t wcycle)
221 real e_q, e_lj, dvdl_q, dvdl_lj;
222 float cycles_ppdpme, cycles_seppme;
224 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
225 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
227 /* In case of node-splitting, the PP nodes receive the long-range
228 * forces, virial and energy from the PME nodes here.
230 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
233 gmx_pme_receive_f(fr->pmePpCommGpu.get(), cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
234 useGpuPmePpComms, receivePmeForceToGpu, &cycles_seppme);
235 enerd->term[F_COUL_RECIP] += e_q;
236 enerd->term[F_LJ_RECIP] += e_lj;
237 enerd->dvdl_lin[efptCOUL] += dvdl_q;
238 enerd->dvdl_lin[efptVDW] += dvdl_lj;
242 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
244 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
247 static void print_large_forces(FILE* fp,
252 ArrayRef<const RVec> x,
253 ArrayRef<const RVec> f)
255 real force2Tolerance = gmx::square(forceTolerance);
256 gmx::index numNonFinite = 0;
257 for (int i = 0; i < md->homenr; i++)
259 real force2 = norm2(f[i]);
260 bool nonFinite = !std::isfinite(force2);
261 if (force2 >= force2Tolerance || nonFinite)
263 fprintf(fp, "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n", step,
264 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
271 if (numNonFinite > 0)
273 /* Note that with MPI this fatal call on one rank might interrupt
274 * the printing on other ranks. But we can only avoid that with
275 * an expensive MPI barrier that we would need at each step.
277 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
281 //! When necessary, spreads forces on vsites and computes the virial for \p forceOutputs->forceWithShiftForces()
282 static void postProcessForceWithShiftForces(t_nrnb* nrnb,
283 gmx_wallcycle_t wcycle,
285 ArrayRef<const RVec> x,
286 ForceOutputs* forceOutputs,
288 const t_mdatoms& mdatoms,
289 const t_forcerec& fr,
290 gmx::VirtualSitesHandler* vsite,
291 const StepWorkload& stepWork)
293 ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
295 /* If we have NoVirSum forces, but we do not calculate the virial,
296 * we later sum the forceWithShiftForces buffer together with
297 * the noVirSum buffer and spread the combined vsite forces at once.
299 if (vsite && (!forceOutputs->haveForceWithVirial() || stepWork.computeVirial))
301 using VirialHandling = gmx::VirtualSitesHandler::VirialHandling;
303 auto f = forceWithShiftForces.force();
304 auto fshift = forceWithShiftForces.shiftForces();
305 const VirialHandling virialHandling =
306 (stepWork.computeVirial ? VirialHandling::Pbc : VirialHandling::None);
307 vsite->spreadForces(x, f, virialHandling, fshift, nullptr, nrnb, box, wcycle);
308 forceWithShiftForces.haveSpreadVsiteForces() = true;
311 if (stepWork.computeVirial)
313 /* Calculation of the virial must be done after vsites! */
314 calc_virial(0, mdatoms.homenr, as_rvec_array(x.data()), forceWithShiftForces, vir_force,
315 box, nrnb, &fr, fr.pbcType);
319 //! Spread, compute virial for and sum forces, when necessary
320 static void postProcessForces(const t_commrec* cr,
323 gmx_wallcycle_t wcycle,
325 ArrayRef<const RVec> x,
326 ForceOutputs* forceOutputs,
328 const t_mdatoms* mdatoms,
329 const t_forcerec* fr,
330 gmx::VirtualSitesHandler* vsite,
331 const StepWorkload& stepWork)
333 // Extract the final output force buffer, which is also the buffer for forces with shift forces
334 ArrayRef<RVec> f = forceOutputs->forceWithShiftForces().force();
336 if (forceOutputs->haveForceWithVirial())
338 auto& forceWithVirial = forceOutputs->forceWithVirial();
342 /* Spread the mesh force on virtual sites to the other particles...
343 * This is parallellized. MPI communication is performed
344 * if the constructing atoms aren't local.
346 GMX_ASSERT(!stepWork.computeVirial || f.data() != forceWithVirial.force_.data(),
347 "We need separate force buffers for shift and virial forces when "
348 "computing the virial");
349 GMX_ASSERT(!stepWork.computeVirial
350 || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
351 "We should spread the force with shift forces separately when computing "
353 const gmx::VirtualSitesHandler::VirialHandling virialHandling =
354 (stepWork.computeVirial ? gmx::VirtualSitesHandler::VirialHandling::NonLinear
355 : gmx::VirtualSitesHandler::VirialHandling::None);
356 matrix virial = { { 0 } };
357 vsite->spreadForces(x, forceWithVirial.force_, virialHandling, {}, virial, nrnb, box, wcycle);
358 forceWithVirial.addVirialContribution(virial);
361 if (stepWork.computeVirial)
363 /* Now add the forces, this is local */
364 sum_forces(f, forceWithVirial.force_);
366 /* Add the direct virial contributions */
368 forceWithVirial.computeVirial_,
369 "forceWithVirial should request virial computation when we request the virial");
370 m_add(vir_force, forceWithVirial.getVirial(), vir_force);
374 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
380 GMX_ASSERT(vsite == nullptr || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
381 "We should have spread the vsite forces (earlier)");
384 if (fr->print_force >= 0)
386 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
390 static void do_nb_verlet(t_forcerec* fr,
391 const interaction_const_t* ic,
392 gmx_enerdata_t* enerd,
393 const StepWorkload& stepWork,
394 const InteractionLocality ilocality,
398 gmx_wallcycle_t wcycle)
400 if (!stepWork.computeNonbondedForces)
402 /* skip non-bonded calculation */
406 nonbonded_verlet_t* nbv = fr->nbv.get();
408 /* GPU kernel launch overhead is already timed separately */
411 /* When dynamic pair-list pruning is requested, we need to prune
412 * at nstlistPrune steps.
414 if (nbv->isDynamicPruningStepCpu(step))
416 /* Prune the pair-list beyond fr->ic->rlistPrune using
417 * the current coordinates of the atoms.
419 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
420 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
421 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
425 nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
428 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
430 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, v.ssize());
432 /* Note that we would like to avoid this conditional by putting it
433 * into the omp pragma instead, but then we still take the full
434 * omp parallel for overhead (at least with gcc5).
436 if (!useOpenmpThreading || nth == 1)
445 #pragma omp parallel for num_threads(nth) schedule(static)
446 for (gmx::index i = 0; i < v.ssize(); i++)
453 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
455 * \param groupOptions Group options, containing T-coupling options
457 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
459 real nrdfCoupled = 0;
460 real nrdfUncoupled = 0;
461 real kineticEnergy = 0;
462 for (int g = 0; g < groupOptions.ngtc; g++)
464 if (groupOptions.tau_t[g] >= 0)
466 nrdfCoupled += groupOptions.nrdf[g];
467 kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * BOLTZ;
471 nrdfUncoupled += groupOptions.nrdf[g];
475 /* This conditional with > also catches nrdf=0 */
476 if (nrdfCoupled > nrdfUncoupled)
478 return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
486 /*! \brief This routine checks that the potential energy is finite.
488 * Always checks that the potential energy is finite. If step equals
489 * inputrec.init_step also checks that the magnitude of the potential energy
490 * is reasonable. Terminates with a fatal error when a check fails.
491 * Note that passing this check does not guarantee finite forces,
492 * since those use slightly different arithmetics. But in most cases
493 * there is just a narrow coordinate range where forces are not finite
494 * and energies are finite.
496 * \param[in] step The step number, used for checking and printing
497 * \param[in] enerd The energy data; the non-bonded group energies need to be added to
498 * enerd.term[F_EPOT] before calling this routine \param[in] inputrec The input record
500 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
502 /* Threshold valid for comparing absolute potential energy against
503 * the kinetic energy. Normally one should not consider absolute
504 * potential energy values, but with a factor of one million
505 * we should never get false positives.
507 constexpr real c_thresholdFactor = 1e6;
509 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
510 real averageKineticEnergy = 0;
511 /* We only check for large potential energy at the initial step,
512 * because that is by far the most likely step for this too occur
513 * and because computing the average kinetic energy is not free.
514 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
515 * before they become NaN.
517 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
519 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
522 if (energyIsNotFinite
523 || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
528 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
529 "contributions to the energy are %g and %g, respectively. A %s potential energy "
530 "can be caused by overlapping interactions in bonded interactions or very large%s "
531 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
532 "configuration, incorrect interactions or parameters in the topology.",
533 step, enerd.term[F_EPOT], energyIsNotFinite ? "not finite" : "extremely high",
534 enerd.term[F_LJ], enerd.term[F_COUL_SR],
535 energyIsNotFinite ? "non-finite" : "very high", energyIsNotFinite ? " or Nan" : "");
539 /*! \brief Return true if there are special forces computed this step.
541 * The conditionals exactly correspond to those in computeSpecialForces().
543 static bool haveSpecialForces(const t_inputrec& inputrec,
544 const gmx::ForceProviders& forceProviders,
545 const pull_t* pull_work,
546 const bool computeForces,
550 return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
551 (inputrec.bPull && pull_have_potential(pull_work)) || // pull
552 inputrec.bRot || // enforced rotation
553 (ed != nullptr) || // flooding
554 (inputrec.bIMD && computeForces)); // IMD
557 /*! \brief Compute forces and/or energies for special algorithms
559 * The intention is to collect all calls to algorithms that compute
560 * forces on local atoms only and that do not contribute to the local
561 * virial sum (but add their virial contribution separately).
562 * Eventually these should likely all become ForceProviders.
563 * Within this function the intention is to have algorithms that do
564 * global communication at the end, so global barriers within the MD loop
565 * are as close together as possible.
567 * \param[in] fplog The log file
568 * \param[in] cr The communication record
569 * \param[in] inputrec The input record
570 * \param[in] awh The Awh module (nullptr if none in use).
571 * \param[in] enforcedRotation Enforced rotation module.
572 * \param[in] imdSession The IMD session
573 * \param[in] pull_work The pull work structure.
574 * \param[in] step The current MD step
575 * \param[in] t The current time
576 * \param[in,out] wcycle Wallcycle accounting struct
577 * \param[in,out] forceProviders Pointer to a list of force providers
578 * \param[in] box The unit cell
579 * \param[in] x The coordinates
580 * \param[in] mdatoms Per atom properties
581 * \param[in] lambda Array of free-energy lambda values
582 * \param[in] stepWork Step schedule flags
583 * \param[in,out] forceWithVirialMtsLevel0 Force and virial for MTS level0 forces
584 * \param[in,out] forceWithVirialMtsLevel1 Force and virial for MTS level1 forces, can be nullptr
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* forceWithVirialMtsLevel0,
609 gmx::ForceWithVirial* forceWithVirialMtsLevel1,
610 gmx_enerdata_t* enerd,
612 bool didNeighborSearch)
614 /* NOTE: Currently all ForceProviders only provide forces.
615 * When they also provide energies, remove this conditional.
617 if (stepWork.computeForces)
619 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
620 gmx::ForceProviderOutput forceProviderOutput(forceWithVirialMtsLevel0, enerd);
622 /* Collect forces from modules */
623 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
626 if (inputrec->bPull && pull_have_potential(pull_work))
628 const int mtsLevel = forceGroupMtsLevel(inputrec->mtsLevels, gmx::MtsForceGroups::Pull);
629 if (mtsLevel == 0 || stepWork.computeSlowForces)
631 auto& forceWithVirial = (mtsLevel == 0) ? forceWithVirialMtsLevel0 : forceWithVirialMtsLevel1;
632 pull_potential_wrapper(cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work,
633 lambda.data(), t, wcycle);
638 const int mtsLevel = forceGroupMtsLevel(inputrec->mtsLevels, gmx::MtsForceGroups::Pull);
639 if (mtsLevel == 0 || stepWork.computeSlowForces)
641 const bool needForeignEnergyDifferences = awh->needForeignEnergyDifferences(step);
642 std::vector<double> foreignLambdaDeltaH, foreignLambdaDhDl;
643 if (needForeignEnergyDifferences)
645 enerd->foreignLambdaTerms.finalizePotentialContributions(enerd->dvdl_lin, lambda,
647 std::tie(foreignLambdaDeltaH, foreignLambdaDhDl) = enerd->foreignLambdaTerms.getTerms(cr);
650 auto& forceWithVirial = (mtsLevel == 0) ? forceWithVirialMtsLevel0 : forceWithVirialMtsLevel1;
651 enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(
652 inputrec->pbcType, mdatoms->massT, foreignLambdaDeltaH, foreignLambdaDhDl, box,
653 forceWithVirial, t, step, wcycle, fplog);
657 rvec* f = as_rvec_array(forceWithVirialMtsLevel0->force_.data());
659 /* Add the forces from enforced rotation potentials (if any) */
662 wallcycle_start(wcycle, ewcROTadd);
663 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
664 wallcycle_stop(wcycle, ewcROTadd);
669 /* Note that since init_edsam() is called after the initialization
670 * of forcerec, edsam doesn't request the noVirSum force buffer.
671 * Thus if no other algorithm (e.g. PME) requires it, the forces
672 * here will contribute to the virial.
674 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
677 /* Add forces from interactive molecular dynamics (IMD), if any */
678 if (inputrec->bIMD && stepWork.computeForces)
680 imdSession->applyForces(f);
684 /*! \brief Launch the prepare_step and spread stages of PME GPU.
686 * \param[in] pmedata The PME structure
687 * \param[in] box The box matrix
688 * \param[in] stepWork Step schedule flags
689 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in the device memory.
690 * \param[in] lambdaQ The Coulomb lambda of the current state.
691 * \param[in] wcycle The wallcycle structure
693 static inline void launchPmeGpuSpread(gmx_pme_t* pmedata,
695 const StepWorkload& stepWork,
696 GpuEventSynchronizer* xReadyOnDevice,
698 gmx_wallcycle_t wcycle)
700 pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
701 pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle, lambdaQ);
704 /*! \brief Launch the FFT and gather stages of PME GPU
706 * This function only implements setting the output forces (no accumulation).
708 * \param[in] pmedata The PME structure
709 * \param[in] lambdaQ The Coulomb lambda of the current system state.
710 * \param[in] wcycle The wallcycle structure
711 * \param[in] stepWork Step schedule flags
713 static void launchPmeGpuFftAndGather(gmx_pme_t* pmedata,
715 gmx_wallcycle_t wcycle,
716 const gmx::StepWorkload& stepWork)
718 pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
719 pme_gpu_launch_gather(pmedata, wcycle, lambdaQ);
723 * Polling wait for either of the PME or nonbonded GPU tasks.
725 * Instead of a static order in waiting for GPU tasks, this function
726 * polls checking which of the two tasks completes first, and does the
727 * associated force buffer reduction overlapped with the other task.
728 * By doing that, unlike static scheduling order, it can always overlap
729 * one of the reductions, regardless of the GPU task completion order.
731 * \param[in] nbv Nonbonded verlet structure
732 * \param[in,out] pmedata PME module data
733 * \param[in,out] forceOutputsNonbonded Force outputs for the non-bonded forces and shift forces
734 * \param[in,out] forceOutputsPme Force outputs for the PME forces and virial
735 * \param[in,out] enerd Energy data structure results are reduced into
736 * \param[in] lambdaQ The Coulomb lambda of the current system state.
737 * \param[in] stepWork Step schedule flags
738 * \param[in] wcycle The wallcycle structure
740 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
742 gmx::ForceOutputs* forceOutputsNonbonded,
743 gmx::ForceOutputs* forceOutputsPme,
744 gmx_enerdata_t* enerd,
746 const StepWorkload& stepWork,
747 gmx_wallcycle_t wcycle)
749 bool isPmeGpuDone = false;
750 bool isNbGpuDone = false;
752 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
754 while (!isPmeGpuDone || !isNbGpuDone)
758 GpuTaskCompletion completionType =
759 (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
760 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, stepWork, wcycle,
761 &forceOutputsPme->forceWithVirial(), enerd,
762 lambdaQ, completionType);
767 auto& forceBuffersNonbonded = forceOutputsNonbonded->forceWithShiftForces();
768 GpuTaskCompletion completionType =
769 (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
770 isNbGpuDone = Nbnxm::gpu_try_finish_task(
771 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
772 enerd->grpp.ener[egCOULSR].data(), forceBuffersNonbonded.shiftForces(),
773 completionType, wcycle);
777 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceBuffersNonbonded.force());
783 /*! \brief Set up the different force buffers; also does clearing.
785 * \param[in] forceHelperBuffers Helper force buffers
786 * \param[in] force force array
787 * \param[in] stepWork Step schedule flags
788 * \param[out] wcycle wallcycle recording structure
790 * \returns Cleared force output structure
792 static ForceOutputs setupForceOutputs(ForceHelperBuffers* forceHelperBuffers,
793 gmx::ArrayRefWithPadding<gmx::RVec> force,
794 const StepWorkload& stepWork,
795 gmx_wallcycle_t wcycle)
797 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
799 /* NOTE: We assume fr->shiftForces is all zeros here */
800 gmx::ForceWithShiftForces forceWithShiftForces(force, stepWork.computeVirial,
801 forceHelperBuffers->shiftForces());
803 if (stepWork.computeForces)
805 /* Clear the short- and long-range forces */
806 clearRVecs(forceWithShiftForces.force(), true);
808 /* Clear the shift forces */
809 clearRVecs(forceWithShiftForces.shiftForces(), false);
812 /* If we need to compute the virial, we might need a separate
813 * force buffer for algorithms for which the virial is calculated
814 * directly, such as PME. Otherwise, forceWithVirial uses the
815 * the same force (f in legacy calls) buffer as other algorithms.
817 const bool useSeparateForceWithVirialBuffer =
818 (stepWork.computeForces
819 && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
820 /* forceWithVirial uses the local atom range only */
821 gmx::ForceWithVirial forceWithVirial(
822 useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
823 : force.unpaddedArrayRef(),
824 stepWork.computeVirial);
826 if (useSeparateForceWithVirialBuffer)
828 /* TODO: update comment
829 * We only compute forces on local atoms. Note that vsites can
830 * spread to non-local atoms, but that part of the buffer is
831 * cleared separately in the vsite spreading code.
833 clearRVecs(forceWithVirial.force_, true);
836 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
838 return ForceOutputs(forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(),
843 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
845 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec& inputrec,
846 const t_forcerec& fr,
847 const pull_t* pull_work,
849 const t_mdatoms& mdatoms,
850 const SimulationWorkload& simulationWork,
851 const StepWorkload& stepWork)
853 DomainLifetimeWorkload domainWork;
854 // Note that haveSpecialForces is constant over the whole run
855 domainWork.haveSpecialForces =
856 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
857 domainWork.haveCpuListedForceWork = false;
858 domainWork.haveCpuBondedWork = false;
859 for (const auto& listedForces : fr.listedForces)
861 if (listedForces.haveCpuListedForces(*fr.fcdata))
863 domainWork.haveCpuListedForceWork = true;
865 if (listedForces.haveCpuBondeds())
867 domainWork.haveCpuBondedWork = true;
870 domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
871 // Note that haveFreeEnergyWork is constant over the whole run
872 domainWork.haveFreeEnergyWork = (fr.efep != efepNO && mdatoms.nPerturbed != 0);
873 // We assume we have local force work if there are CPU
874 // force tasks including PME or nonbondeds.
875 domainWork.haveCpuLocalForceWork =
876 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
877 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
878 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
883 /*! \brief Set up force flag stuct from the force bitmask.
885 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
886 * \param[in] mtsLevels The multiple time-stepping levels, either empty or 2 levels
887 * \param[in] step The current MD step
888 * \param[in] simulationWork Simulation workload description.
889 * \param[in] rankHasPmeDuty If this rank computes PME.
891 * \returns New Stepworkload description.
893 static StepWorkload setupStepWorkload(const int legacyFlags,
894 ArrayRef<const gmx::MtsLevel> mtsLevels,
896 const SimulationWorkload& simulationWork,
897 const bool rankHasPmeDuty)
899 GMX_ASSERT(mtsLevels.empty() || mtsLevels.size() == 2, "Expect 0 or 2 MTS levels");
900 const bool computeSlowForces = (mtsLevels.empty() || step % mtsLevels[1].stepFactor == 0);
903 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
904 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
905 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
906 flags.computeSlowForces = computeSlowForces;
907 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
908 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
909 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
910 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
911 flags.computeNonbondedForces =
912 ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded
913 && !(simulationWork.computeNonbondedAtMtsLevel1 && !computeSlowForces);
914 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
916 if (simulationWork.useGpuBufferOps)
918 GMX_ASSERT(simulationWork.useGpuNonbonded,
919 "Can only offload buffer ops if nonbonded computation is also offloaded");
921 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
922 // on virial steps the CPU reduction path is taken
923 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
924 flags.useGpuPmeFReduction = flags.computeSlowForces && flags.useGpuFBufferOps && simulationWork.useGpuPme
925 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication);
926 flags.useGpuXHalo = simulationWork.useGpuHaloExchange;
927 flags.useGpuFHalo = simulationWork.useGpuHaloExchange && flags.useGpuFBufferOps;
933 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
935 * TODO: eliminate \p useGpuPmeOnThisRank when this is
936 * incorporated in DomainLifetimeWorkload.
938 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
939 gmx::GpuBonded* gpuBonded,
941 gmx_enerdata_t* enerd,
942 const gmx::MdrunScheduleWorkload& runScheduleWork,
943 bool useGpuPmeOnThisRank,
945 gmx_wallcycle_t wcycle)
947 if (runScheduleWork.simulationWork.useGpuNonbonded && runScheduleWork.stepWork.computeNonbondedForces)
949 /* Launch pruning before buffer clearing because the API overhead of the
950 * clear kernel launches can leave the GPU idle while it could be running
953 if (nbv->isDynamicPruningStepGpu(step))
955 nbv->dispatchPruneKernelGpu(step);
958 /* now clear the GPU outputs while we finish the step on the CPU */
959 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
960 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
961 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
962 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
963 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
966 if (useGpuPmeOnThisRank)
968 pme_gpu_reinit_computation(pmedata, wcycle);
971 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
973 // in principle this should be included in the DD balancing region,
974 // but generally it is infrequent so we'll omit it for the sake of
976 gpuBonded->waitAccumulateEnergyTerms(enerd);
978 gpuBonded->clearEnergies();
982 //! \brief Data structure to hold dipole-related data and staging arrays
985 //! Dipole staging for fast summing over MPI
986 gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
987 //! Dipole staging for states A and B (index 0 and 1 resp.)
988 gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
992 static void reduceAndUpdateMuTot(DipoleData* dipoleData,
994 const bool haveFreeEnergy,
995 gmx::ArrayRef<const real> lambda,
997 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1001 gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
1002 ddBalanceRegionHandler.reopenRegionCpu();
1004 for (int i = 0; i < 2; i++)
1006 for (int j = 0; j < DIM; j++)
1008 dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
1012 if (!haveFreeEnergy)
1014 copy_rvec(dipoleData->muStateAB[0], muTotal);
1018 for (int j = 0; j < DIM; j++)
1020 muTotal[j] = (1.0 - lambda[efptCOUL]) * dipoleData->muStateAB[0][j]
1021 + lambda[efptCOUL] * dipoleData->muStateAB[1][j];
1026 /*! \brief Combines MTS level0 and level1 force buffes into a full and MTS-combined force buffer.
1028 * \param[in] numAtoms The number of atoms to combine forces for
1029 * \param[in,out] forceMtsLevel0 Input: F_level0, output: F_level0 + F_level1
1030 * \param[in,out] forceMts Input: F_level1, output: F_level0 + mtsFactor * F_level1
1031 * \param[in] mtsFactor The factor between the level0 and level1 time step
1033 static void combineMtsForces(const int numAtoms,
1034 ArrayRef<RVec> forceMtsLevel0,
1035 ArrayRef<RVec> forceMts,
1036 const real mtsFactor)
1038 const int gmx_unused numThreads = gmx_omp_nthreads_get(emntDefault);
1039 #pragma omp parallel for num_threads(numThreads) schedule(static)
1040 for (int i = 0; i < numAtoms; i++)
1042 const RVec forceMtsLevel0Tmp = forceMtsLevel0[i];
1043 forceMtsLevel0[i] += forceMts[i];
1044 forceMts[i] = forceMtsLevel0Tmp + mtsFactor * forceMts[i];
1048 /*! \brief Setup for the local and non-local GPU force reductions:
1049 * reinitialization plus the registration of forces and dependencies.
1051 * \param [in] runScheduleWork Schedule workload flag structure
1052 * \param [in] cr Communication record object
1053 * \param [in] fr Force record object
1055 static void setupGpuForceReductions(gmx::MdrunScheduleWorkload* runScheduleWork,
1056 const t_commrec* cr,
1060 nonbonded_verlet_t* nbv = fr->nbv.get();
1061 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1063 // (re-)initialize local GPU force reduction
1064 const bool accumulate =
1065 runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr);
1066 const int atomStart = 0;
1067 fr->gpuForceReduction[gmx::AtomLocality::Local]->reinit(
1068 stateGpu->getForces(), nbv->getNumAtoms(AtomLocality::Local), nbv->getGridIndices(),
1069 atomStart, accumulate, stateGpu->fReducedOnDevice());
1071 // register forces and add dependencies
1072 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerNbnxmForce(nbv->getGpuForces());
1074 if (runScheduleWork->simulationWork.useGpuPme
1075 && (thisRankHasDuty(cr, DUTY_PME) || runScheduleWork->simulationWork.useGpuPmePpCommunication))
1077 void* forcePtr = thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1078 : // PME force buffer on same GPU
1079 fr->pmePpCommGpu->getGpuForceStagingPtr(); // buffer received from other GPU
1080 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerRvecForce(forcePtr);
1082 GpuEventSynchronizer* const pmeSynchronizer =
1083 (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1084 : // PME force buffer on same GPU
1085 fr->pmePpCommGpu->getForcesReadySynchronizer()); // buffer received from other GPU
1086 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(pmeSynchronizer);
1089 if ((runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr))
1090 && !runScheduleWork->simulationWork.useGpuHaloExchange)
1092 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1093 stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::Local, true));
1096 if (runScheduleWork->simulationWork.useGpuHaloExchange)
1098 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1099 cr->dd->gpuHaloExchange[0][0]->getForcesReadyOnDeviceEvent());
1102 if (havePPDomainDecomposition(cr))
1104 // (re-)initialize non-local GPU force reduction
1105 const bool accumulate = runScheduleWork->domainWork.haveCpuBondedWork
1106 || runScheduleWork->domainWork.haveFreeEnergyWork;
1107 const int atomStart = dd_numHomeAtoms(*cr->dd);
1108 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->reinit(
1109 stateGpu->getForces(), nbv->getNumAtoms(AtomLocality::NonLocal),
1110 nbv->getGridIndices(), atomStart, accumulate);
1112 // register forces and add dependencies
1113 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->registerNbnxmForce(nbv->getGpuForces());
1114 if (runScheduleWork->domainWork.haveCpuBondedWork || runScheduleWork->domainWork.haveFreeEnergyWork)
1116 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->addDependency(
1117 stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::NonLocal, true));
1123 void do_force(FILE* fplog,
1124 const t_commrec* cr,
1125 const gmx_multisim_t* ms,
1126 const t_inputrec* inputrec,
1128 gmx_enfrot* enforcedRotation,
1129 gmx::ImdSession* imdSession,
1133 gmx_wallcycle_t wcycle,
1134 const gmx_localtop_t* top,
1136 gmx::ArrayRefWithPadding<gmx::RVec> x,
1138 gmx::ForceBuffersView* forceView,
1140 const t_mdatoms* mdatoms,
1141 gmx_enerdata_t* enerd,
1142 gmx::ArrayRef<const real> lambda,
1144 gmx::MdrunScheduleWorkload* runScheduleWork,
1145 gmx::VirtualSitesHandler* vsite,
1150 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1152 auto force = forceView->forceWithPadding();
1153 GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1154 "The size of the force buffer should be at least the number of atoms to compute "
1157 nonbonded_verlet_t* nbv = fr->nbv.get();
1158 interaction_const_t* ic = fr->ic;
1160 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1162 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1164 runScheduleWork->stepWork = setupStepWorkload(legacyFlags, inputrec->mtsLevels, step,
1165 simulationWork, thisRankHasDuty(cr, DUTY_PME));
1166 const StepWorkload& stepWork = runScheduleWork->stepWork;
1168 const bool useGpuPmeOnThisRank =
1169 simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces;
1171 /* At a search step we need to start the first balancing region
1172 * somewhere early inside the step after communication during domain
1173 * decomposition (and not during the previous step as usual).
1175 if (stepWork.doNeighborSearch)
1177 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1180 clear_mat(vir_force);
1182 if (fr->pbcType != PbcType::No)
1184 /* Compute shift vectors every step,
1185 * because of pressure coupling or box deformation!
1187 if (stepWork.haveDynamicBox && stepWork.stateChanged)
1189 calc_shifts(box, fr->shift_vec);
1192 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1193 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1196 put_atoms_in_box_omp(fr->pbcType, box, x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1197 gmx_omp_nthreads_get(emntDefault));
1198 inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1202 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1204 const bool pmeSendCoordinatesFromGpu =
1205 GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1206 const bool reinitGpuPmePpComms =
1207 GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1209 const auto localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1210 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1211 AtomLocality::Local, simulationWork, stepWork)
1214 // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1215 // Otherwise the send will occur after H2D coordinate transfer.
1216 if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu && stepWork.computeSlowForces)
1218 /* Send particle coordinates to the pme nodes */
1219 if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1221 GMX_RELEASE_ASSERT(false,
1222 "GPU update and separate PME ranks are only supported with GPU "
1223 "direct communication!");
1224 // TODO: when this code-path becomes supported add:
1225 // stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1228 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1229 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1230 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1231 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1234 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1235 // The local coordinates can be copied right away.
1236 // NOTE: Consider moving this copy to right after they are updated and constrained,
1237 // if the later is not offloaded.
1238 if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1240 if (stepWork.doNeighborSearch)
1242 // TODO refactor this to do_md, after partitioning.
1243 stateGpu->reinit(mdatoms->homenr,
1244 cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1245 if (useGpuPmeOnThisRank)
1247 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1248 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1251 // We need to copy coordinates when:
1252 // 1. Update is not offloaded
1253 // 2. The buffers were reinitialized on search step
1254 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1256 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1257 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1261 // Copy coordinate from the GPU if update is on the GPU and there
1262 // are forces to be computed on the CPU, or for the computation of
1263 // virial, or if host-side data will be transferred from this task
1264 // to a remote task for halo exchange or PME-PP communication. At
1265 // search steps the current coordinates are already on the host,
1266 // hence copy is not needed.
1267 const bool haveHostPmePpComms =
1268 !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1270 GMX_ASSERT(simulationWork.useGpuHaloExchange
1271 == ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange[0].empty())),
1272 "The GPU halo exchange is active, but it has not been constructed.");
1273 const bool haveHostHaloExchangeComms =
1274 havePPDomainDecomposition(cr) && !simulationWork.useGpuHaloExchange;
1276 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1277 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1278 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1279 || haveHostPmePpComms || haveHostHaloExchangeComms))
1281 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1282 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1283 haveCopiedXFromGpu = true;
1286 // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1287 // Otherwise the send will occur before the H2D coordinate transfer.
1288 if (!thisRankHasDuty(cr, DUTY_PME) && pmeSendCoordinatesFromGpu)
1290 /* Send particle coordinates to the pme nodes */
1291 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1292 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1293 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1294 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1297 if (useGpuPmeOnThisRank)
1299 launchPmeGpuSpread(fr->pmedata, box, stepWork, localXReadyOnDevice, lambda[efptCOUL], wcycle);
1302 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1304 /* do gridding for pair search */
1305 if (stepWork.doNeighborSearch)
1307 if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1309 fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1313 // - vzero is constant, do we need to pass it?
1314 // - box_diag should be passed directly to nbnxn_put_on_grid
1320 box_diag[XX] = box[XX][XX];
1321 box_diag[YY] = box[YY][YY];
1322 box_diag[ZZ] = box[ZZ][ZZ];
1324 wallcycle_start(wcycle, ewcNS);
1325 if (!DOMAINDECOMP(cr))
1327 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1328 nbnxn_put_on_grid(nbv, box, 0, vzero, box_diag, nullptr, { 0, mdatoms->homenr }, -1,
1329 fr->cginfo, x.unpaddedArrayRef(), 0, nullptr);
1330 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1334 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1335 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1336 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1339 nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1340 gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr), fr->cginfo);
1342 wallcycle_stop(wcycle, ewcNS);
1344 /* initialize the GPU nbnxm atom data and bonded data structures */
1345 if (simulationWork.useGpuNonbonded)
1347 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1349 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1350 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1351 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1355 /* Now we put all atoms on the grid, we can assign bonded
1356 * interactions to the GPU, where the grid order is
1357 * needed. Also the xq, f and fshift device buffers have
1358 * been reallocated if needed, so the bonded code can
1359 * learn about them. */
1360 // TODO the xq, f, and fshift buffers are now shared
1361 // resources, so they should be maintained by a
1362 // higher-level object than the nb module.
1363 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(
1364 nbv->getGridIndices(), top->idef, Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1365 Nbnxm::gpu_get_f(nbv->gpu_nbv), Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1367 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1370 // Need to run after the GPU-offload bonded interaction lists
1371 // are set up to be able to determine whether there is bonded work.
1372 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1373 *inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1375 wallcycle_start_nocount(wcycle, ewcNS);
1376 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1377 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1378 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1380 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1382 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1383 wallcycle_stop(wcycle, ewcNS);
1385 if (stepWork.useGpuXBufferOps)
1387 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1390 if (simulationWork.useGpuBufferOps)
1392 setupGpuForceReductions(runScheduleWork, cr, fr);
1395 else if (!EI_TPI(inputrec->eI) && stepWork.computeNonbondedForces)
1397 if (stepWork.useGpuXBufferOps)
1399 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1400 nbv->convertCoordinatesGpu(AtomLocality::Local, false, stateGpu->getCoordinates(),
1401 localXReadyOnDevice);
1405 if (simulationWork.useGpuUpdate)
1407 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1408 GMX_ASSERT(haveCopiedXFromGpu,
1409 "a wait should only be triggered if copy has been scheduled");
1410 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1412 nbv->convertCoordinates(AtomLocality::Local, false, x.unpaddedArrayRef());
1416 if (simulationWork.useGpuNonbonded && (stepWork.computeNonbondedForces || domainWork.haveGpuBondedWork))
1418 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1420 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1422 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1423 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1424 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1426 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1428 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1429 // with X buffer ops offloaded to the GPU on all but the search steps
1431 // bonded work not split into separate local and non-local, so with DD
1432 // we can only launch the kernel after non-local coordinates have been received.
1433 if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1435 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1436 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1437 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1440 /* launch local nonbonded work on GPU */
1441 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1442 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1443 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1444 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1447 if (useGpuPmeOnThisRank)
1449 // In PME GPU and mixed mode we launch FFT / gather after the
1450 // X copy/transform to allow overlap as well as after the GPU NB
1451 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1452 // the nonbonded kernel.
1453 launchPmeGpuFftAndGather(fr->pmedata, lambda[efptCOUL], wcycle, stepWork);
1456 /* Communicate coordinates and sum dipole if necessary +
1457 do non-local pair search */
1458 if (havePPDomainDecomposition(cr))
1460 if (stepWork.doNeighborSearch)
1462 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1463 wallcycle_start_nocount(wcycle, ewcNS);
1464 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1465 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1466 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1468 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1469 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1470 wallcycle_stop(wcycle, ewcNS);
1471 // TODO refactor this GPU halo exchange re-initialisation
1472 // to location in do_md where GPU halo exchange is
1473 // constructed at partitioning, after above stateGpu
1474 // re-initialization has similarly been refactored
1475 if (simulationWork.useGpuHaloExchange)
1477 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1482 if (stepWork.useGpuXHalo)
1484 // The following must be called after local setCoordinates (which records an event
1485 // when the coordinate data has been copied to the device).
1486 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1488 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1490 // non-local part of coordinate buffer must be copied back to host for CPU work
1491 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1496 // Note: GPU update + DD without direct communication is not supported,
1497 // a waitCoordinatesReadyOnHost() should be issued if it will be.
1498 GMX_ASSERT(!simulationWork.useGpuUpdate,
1499 "GPU update is not supported with CPU halo exchange");
1500 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1503 if (stepWork.useGpuXBufferOps)
1505 if (!useGpuPmeOnThisRank && !stepWork.useGpuXHalo)
1507 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1509 nbv->convertCoordinatesGpu(AtomLocality::NonLocal, false, stateGpu->getCoordinates(),
1510 stateGpu->getCoordinatesReadyOnDeviceEvent(
1511 AtomLocality::NonLocal, simulationWork, stepWork));
1515 nbv->convertCoordinates(AtomLocality::NonLocal, false, x.unpaddedArrayRef());
1519 if (simulationWork.useGpuNonbonded)
1521 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1523 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1525 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1526 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1527 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1530 if (domainWork.haveGpuBondedWork)
1532 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1533 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1534 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1537 /* launch non-local nonbonded tasks on GPU */
1538 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1539 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1541 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1543 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1547 if (simulationWork.useGpuNonbonded && stepWork.computeNonbondedForces)
1549 /* launch D2H copy-back F */
1550 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1551 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1553 if (havePPDomainDecomposition(cr))
1555 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1557 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1558 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1560 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1562 fr->gpuBonded->launchEnergyTransfer();
1564 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1567 gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1568 if (fr->wholeMoleculeTransform)
1570 xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1573 DipoleData dipoleData;
1575 if (simulationWork.computeMuTot)
1577 const int start = 0;
1579 /* Calculate total (local) dipole moment in a temporary common array.
1580 * This makes it possible to sum them over nodes faster.
1582 gmx::ArrayRef<const gmx::RVec> xRef =
1583 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1584 calc_mu(start, mdatoms->homenr, xRef, mdatoms->chargeA, mdatoms->chargeB,
1585 mdatoms->nChargePerturbed, dipoleData.muStaging[0], dipoleData.muStaging[1]);
1587 reduceAndUpdateMuTot(&dipoleData, cr, (fr->efep != efepNO), lambda, muTotal, ddBalanceRegionHandler);
1590 /* Reset energies */
1591 reset_enerdata(enerd);
1593 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1595 wallcycle_start(wcycle, ewcPPDURINGPME);
1596 dd_force_flop_start(cr->dd, nrnb);
1599 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1600 // this wait ensures that the D2H transfer is complete.
1601 if ((simulationWork.useGpuUpdate)
1602 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1604 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1609 wallcycle_start(wcycle, ewcROT);
1610 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step,
1611 stepWork.doNeighborSearch);
1612 wallcycle_stop(wcycle, ewcROT);
1615 /* Start the force cycle counter.
1616 * Note that a different counter is used for dynamic load balancing.
1618 wallcycle_start(wcycle, ewcFORCE);
1620 /* Set up and clear force outputs:
1621 * forceOutMtsLevel0: everything except what is in the other two outputs
1622 * forceOutMtsLevel1: PME-mesh and listed-forces group 1
1623 * forceOutNonbonded: non-bonded forces
1624 * Without multiple time stepping all point to the same object.
1625 * With multiple time-stepping the use is different for MTS fast (level0 only) and slow steps.
1627 ForceOutputs forceOutMtsLevel0 =
1628 setupForceOutputs(&fr->forceHelperBuffers[0], force, stepWork, wcycle);
1630 // Force output for MTS combined forces, only set at level1 MTS steps
1631 std::optional<ForceOutputs> forceOutMts =
1632 (fr->useMts && stepWork.computeSlowForces)
1633 ? std::optional(setupForceOutputs(&fr->forceHelperBuffers[1],
1634 forceView->forceMtsCombinedWithPadding(),
1638 ForceOutputs* forceOutMtsLevel1 =
1639 fr->useMts ? (stepWork.computeSlowForces ? &forceOutMts.value() : nullptr) : &forceOutMtsLevel0;
1641 const bool nonbondedAtMtsLevel1 = runScheduleWork->simulationWork.computeNonbondedAtMtsLevel1;
1643 ForceOutputs* forceOutNonbonded = nonbondedAtMtsLevel1 ? forceOutMtsLevel1 : &forceOutMtsLevel0;
1645 if (inputrec->bPull && pull_have_constraint(pull_work))
1647 clear_pull_forces(pull_work);
1650 /* We calculate the non-bonded forces, when done on the CPU, here.
1651 * We do this before calling do_force_lowlevel, because in that
1652 * function, the listed forces are calculated before PME, which
1653 * does communication. With this order, non-bonded and listed
1654 * force calculation imbalance can be balanced out by the domain
1655 * decomposition load balancing.
1658 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1660 if (!useOrEmulateGpuNb)
1662 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1665 if (fr->efep != efepNO && stepWork.computeNonbondedForces)
1667 /* Calculate the local and non-local free energy interactions here.
1668 * Happens here on the CPU both with and without GPU.
1670 nbv->dispatchFreeEnergyKernel(InteractionLocality::Local, fr,
1671 as_rvec_array(x.unpaddedArrayRef().data()),
1672 &forceOutNonbonded->forceWithShiftForces(), *mdatoms,
1673 inputrec->fepvals, lambda, enerd, stepWork, nrnb);
1675 if (havePPDomainDecomposition(cr))
1677 nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal, fr,
1678 as_rvec_array(x.unpaddedArrayRef().data()),
1679 &forceOutNonbonded->forceWithShiftForces(), *mdatoms,
1680 inputrec->fepvals, lambda, enerd, stepWork, nrnb);
1684 if (stepWork.computeNonbondedForces && !useOrEmulateGpuNb)
1686 if (havePPDomainDecomposition(cr))
1688 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1692 if (stepWork.computeForces)
1694 /* Add all the non-bonded force to the normal force array.
1695 * This can be split into a local and a non-local part when overlapping
1696 * communication with calculation with domain decomposition.
1698 wallcycle_stop(wcycle, ewcFORCE);
1699 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All,
1700 forceOutNonbonded->forceWithShiftForces().force());
1701 wallcycle_start_nocount(wcycle, ewcFORCE);
1704 /* If there are multiple fshift output buffers we need to reduce them */
1705 if (stepWork.computeVirial)
1707 /* This is not in a subcounter because it takes a
1708 negligible and constant-sized amount of time */
1709 nbnxn_atomdata_add_nbat_fshift_to_fshift(
1710 *nbv->nbat, forceOutNonbonded->forceWithShiftForces().shiftForces());
1714 // TODO Force flags should include haveFreeEnergyWork for this domain
1715 if (stepWork.useGpuXHalo && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1717 /* Wait for non-local coordinate data to be copied from device */
1718 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1721 // Compute wall interactions, when present.
1722 // Note: should be moved to special forces.
1723 if (inputrec->nwall && stepWork.computeNonbondedForces)
1725 /* foreign lambda component for walls */
1726 real dvdl_walls = do_walls(*inputrec, *fr, box, *mdatoms, x.unpaddedConstArrayRef(),
1727 &forceOutMtsLevel0.forceWithVirial(), lambda[efptVDW],
1728 enerd->grpp.ener[egLJSR].data(), nrnb);
1729 enerd->dvdl_lin[efptVDW] += dvdl_walls;
1732 if (stepWork.computeListedForces)
1734 /* Check whether we need to take into account PBC in listed interactions */
1735 bool needMolPbc = false;
1736 for (const auto& listedForces : fr->listedForces)
1738 if (listedForces.haveCpuListedForces(*fr->fcdata))
1740 needMolPbc = fr->bMolPBC;
1748 /* Since all atoms are in the rectangular or triclinic unit-cell,
1749 * only single box vector shifts (2 in x) are required.
1751 set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1754 for (int mtsIndex = 0; mtsIndex < (fr->useMts && stepWork.computeSlowForces ? 2 : 1); mtsIndex++)
1756 ListedForces& listedForces = fr->listedForces[mtsIndex];
1757 ForceOutputs& forceOut = (mtsIndex == 0 ? forceOutMtsLevel0 : *forceOutMtsLevel1);
1758 listedForces.calculate(
1759 wcycle, box, inputrec->fepvals, cr, ms, x, xWholeMolecules, fr->fcdata.get(),
1760 hist, &forceOut, fr, &pbc, enerd, nrnb, lambda.data(), mdatoms,
1761 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr, stepWork);
1765 if (stepWork.computeSlowForces)
1767 calculateLongRangeNonbondeds(fr, inputrec, cr, nrnb, wcycle, mdatoms,
1768 x.unpaddedConstArrayRef(), &forceOutMtsLevel1->forceWithVirial(),
1769 enerd, box, lambda.data(), as_rvec_array(dipoleData.muStateAB),
1770 stepWork, ddBalanceRegionHandler);
1773 wallcycle_stop(wcycle, ewcFORCE);
1775 // VdW dispersion correction, only computed on master rank to avoid double counting
1776 if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1778 // Calculate long range corrections to pressure and energy
1779 const DispersionCorrection::Correction correction =
1780 fr->dispersionCorrection->calculate(box, lambda[efptVDW]);
1782 if (stepWork.computeEnergy)
1784 enerd->term[F_DISPCORR] = correction.energy;
1785 enerd->term[F_DVDL_VDW] += correction.dvdl;
1786 enerd->dvdl_lin[efptVDW] += correction.dvdl;
1788 if (stepWork.computeVirial)
1790 correction.correctVirial(vir_force);
1791 enerd->term[F_PDISPCORR] = correction.pressure;
1795 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation, imdSession, pull_work, step, t,
1796 wcycle, fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda,
1797 stepWork, &forceOutMtsLevel0.forceWithVirial(),
1798 forceOutMtsLevel1 ? &forceOutMtsLevel1->forceWithVirial() : nullptr, enerd,
1799 ed, stepWork.doNeighborSearch);
1801 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
1802 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
1803 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFHalo),
1804 "The schedule below does not allow for nonbonded MTS with GPU halo exchange");
1805 // Will store the amount of cycles spent waiting for the GPU that
1806 // will be later used in the DLB accounting.
1807 float cycles_wait_gpu = 0;
1808 if (useOrEmulateGpuNb && stepWork.computeNonbondedForces)
1810 auto& forceWithShiftForces = forceOutNonbonded->forceWithShiftForces();
1812 /* wait for non-local forces (or calculate in emulation mode) */
1813 if (havePPDomainDecomposition(cr))
1815 if (simulationWork.useGpuNonbonded)
1817 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1818 nbv->gpu_nbv, stepWork, AtomLocality::NonLocal, enerd->grpp.ener[egLJSR].data(),
1819 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(), wcycle);
1823 wallcycle_start_nocount(wcycle, ewcFORCE);
1824 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes,
1825 step, nrnb, wcycle);
1826 wallcycle_stop(wcycle, ewcFORCE);
1829 if (stepWork.useGpuFBufferOps)
1831 // TODO: move this into DomainLifetimeWorkload, including the second part of the
1832 // condition The bonded and free energy CPU tasks can have non-local force
1833 // contributions which are a dependency for the GPU force reduction.
1834 bool haveNonLocalForceContribInCpuBuffer =
1835 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1837 if (haveNonLocalForceContribInCpuBuffer)
1839 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
1840 AtomLocality::NonLocal);
1844 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->execute();
1846 if (!stepWork.useGpuFHalo)
1848 // copy from GPU input for dd_move_f()
1849 stateGpu->copyForcesFromGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
1850 AtomLocality::NonLocal);
1855 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
1858 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
1860 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
1865 /* Combining the forces for multiple time stepping before the halo exchange, when possible,
1866 * avoids an extra halo exchange (when DD is used) and post-processing step.
1868 const bool combineMtsForcesBeforeHaloExchange =
1869 (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
1870 && (legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0
1871 && !(stepWork.computeVirial || simulationWork.useGpuNonbonded || useGpuPmeOnThisRank));
1872 if (combineMtsForcesBeforeHaloExchange)
1874 const int numAtoms = havePPDomainDecomposition(cr) ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr;
1875 combineMtsForces(numAtoms, force.unpaddedArrayRef(), forceView->forceMtsCombined(),
1876 inputrec->mtsLevels[1].stepFactor);
1879 if (havePPDomainDecomposition(cr))
1881 /* We are done with the CPU compute.
1882 * We will now communicate the non-local forces.
1883 * If we use a GPU this will overlap with GPU work, so in that case
1884 * we do not close the DD force balancing region here.
1886 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1888 if (stepWork.computeForces)
1891 if (stepWork.useGpuFHalo)
1893 if (domainWork.haveCpuLocalForceWork)
1895 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
1896 AtomLocality::Local);
1898 communicateGpuHaloForces(*cr, domainWork.haveCpuLocalForceWork);
1902 if (stepWork.useGpuFBufferOps)
1904 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
1907 // Without MTS or with MTS at slow steps with uncombined forces we need to
1908 // communicate the fast forces
1909 if (!fr->useMts || !combineMtsForcesBeforeHaloExchange)
1911 dd_move_f(cr->dd, &forceOutMtsLevel0.forceWithShiftForces(), wcycle);
1913 // With MTS we need to communicate the slow or combined (in forceOutMtsLevel1) forces
1914 if (fr->useMts && stepWork.computeSlowForces)
1916 dd_move_f(cr->dd, &forceOutMtsLevel1->forceWithShiftForces(), wcycle);
1922 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1923 // an alternating wait/reduction scheme.
1924 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
1925 && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
1926 if (alternateGpuWait)
1928 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, forceOutNonbonded,
1929 forceOutMtsLevel1, enerd, lambda[efptCOUL], stepWork, wcycle);
1932 if (!alternateGpuWait && useGpuPmeOnThisRank)
1934 pme_gpu_wait_and_reduce(fr->pmedata, stepWork, wcycle,
1935 &forceOutMtsLevel1->forceWithVirial(), enerd, lambda[efptCOUL]);
1938 /* Wait for local GPU NB outputs on the non-alternating wait path */
1939 if (!alternateGpuWait && stepWork.computeNonbondedForces && simulationWork.useGpuNonbonded)
1941 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1942 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1943 * but even with a step of 0.1 ms the difference is less than 1%
1946 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
1947 const float waitCycles = Nbnxm::gpu_wait_finish_task(
1948 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
1949 enerd->grpp.ener[egCOULSR].data(),
1950 forceOutNonbonded->forceWithShiftForces().shiftForces(), wcycle);
1952 if (ddBalanceRegionHandler.useBalancingRegion())
1954 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1955 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
1957 /* We measured few cycles, it could be that the kernel
1958 * and transfer finished earlier and there was no actual
1959 * wait time, only API call overhead.
1960 * Then the actual time could be anywhere between 0 and
1961 * cycles_wait_est. We will use half of cycles_wait_est.
1963 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1965 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1969 if (fr->nbv->emulateGpu())
1971 // NOTE: emulation kernel is not included in the balancing region,
1972 // but emulation mode does not target performance anyway
1973 wallcycle_start_nocount(wcycle, ewcFORCE);
1974 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local,
1975 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes, step, nrnb, wcycle);
1976 wallcycle_stop(wcycle, ewcFORCE);
1979 // If on GPU PME-PP comms or GPU update path, receive forces from PME before GPU buffer ops
1980 // TODO refactor this and unify with below default-path call to the same function
1981 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces
1982 && (simulationWork.useGpuPmePpCommunication || simulationWork.useGpuUpdate))
1984 /* In case of node-splitting, the PP nodes receive the long-range
1985 * forces, virial and energy from the PME nodes here.
1987 pme_receive_force_ener(fr, cr, &forceOutMtsLevel1->forceWithVirial(), enerd,
1988 simulationWork.useGpuPmePpCommunication,
1989 stepWork.useGpuPmeFReduction, wcycle);
1993 /* Do the nonbonded GPU (or emulation) force buffer reduction
1994 * on the non-alternating path. */
1995 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
1996 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
1997 if (useOrEmulateGpuNb && !alternateGpuWait)
1999 if (stepWork.useGpuFBufferOps)
2001 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2003 // Flag to specify whether the CPU force buffer has contributions to
2004 // local atoms. This depends on whether there are CPU-based force tasks
2005 // or when DD is active the halo exchange has resulted in contributions
2006 // from the non-local part.
2007 const bool haveLocalForceContribInCpuBuffer =
2008 (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
2010 // TODO: move these steps as early as possible:
2011 // - CPU f H2D should be as soon as all CPU-side forces are done
2012 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
2013 // before the next CPU task that consumes the forces: vsite spread or update)
2014 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
2015 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
2016 // These should be unified.
2017 if (haveLocalForceContribInCpuBuffer && !stepWork.useGpuFHalo)
2019 // Note: AtomLocality::All is used for the non-DD case because, as in this
2020 // case copyForcesToGpu() uses a separate stream, it allows overlap of
2021 // CPU force H2D with GPU force tasks on all streams including those in the
2022 // local stream which would otherwise be implicit dependencies for the
2023 // transfer and would not overlap.
2024 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
2026 stateGpu->copyForcesToGpu(forceWithShift, locality);
2029 if (stepWork.computeNonbondedForces)
2031 fr->gpuForceReduction[gmx::AtomLocality::Local]->execute();
2034 // Copy forces to host if they are needed for update or if virtual sites are enabled.
2035 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
2036 // TODO: When the output flags will be included in step workload, this copy can be combined with the
2037 // copy call done in sim_utils(...) for the output.
2038 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
2039 // they should not be copied in do_md(...) for the output.
2040 if (!simulationWork.useGpuUpdate || vsite)
2042 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
2043 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
2046 else if (stepWork.computeNonbondedForces)
2048 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2049 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
2053 launchGpuEndOfStepTasks(nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork,
2054 useGpuPmeOnThisRank, step, wcycle);
2056 if (DOMAINDECOMP(cr))
2058 dd_force_flop_stop(cr->dd, nrnb);
2061 const bool haveCombinedMtsForces = (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2062 && combineMtsForcesBeforeHaloExchange);
2063 if (stepWork.computeForces)
2065 postProcessForceWithShiftForces(nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutMtsLevel0,
2066 vir_force, *mdatoms, *fr, vsite, stepWork);
2068 if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2070 postProcessForceWithShiftForces(nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1,
2071 vir_force, *mdatoms, *fr, vsite, stepWork);
2075 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
2076 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
2077 && !simulationWork.useGpuUpdate && stepWork.computeSlowForces)
2079 /* In case of node-splitting, the PP nodes receive the long-range
2080 * forces, virial and energy from the PME nodes here.
2082 pme_receive_force_ener(fr, cr, &forceOutMtsLevel1->forceWithVirial(), enerd,
2083 simulationWork.useGpuPmePpCommunication, false, wcycle);
2086 if (stepWork.computeForces)
2088 /* If we don't use MTS or if we already combined the MTS forces before, we only
2089 * need to post-process one ForceOutputs object here, called forceOutCombined,
2090 * otherwise we have to post-process two outputs and then combine them.
2092 ForceOutputs& forceOutCombined = (haveCombinedMtsForces ? forceOutMts.value() : forceOutMtsLevel0);
2093 postProcessForces(cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutCombined,
2094 vir_force, mdatoms, fr, vsite, stepWork);
2096 if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2098 postProcessForces(cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1,
2099 vir_force, mdatoms, fr, vsite, stepWork);
2101 combineMtsForces(mdatoms->homenr, force.unpaddedArrayRef(),
2102 forceView->forceMtsCombined(), inputrec->mtsLevels[1].stepFactor);
2106 if (stepWork.computeEnergy)
2108 /* Compute the final potential energy terms */
2109 accumulatePotentialEnergies(enerd, lambda, inputrec->fepvals);
2111 if (!EI_TPI(inputrec->eI))
2113 checkPotentialEnergyValidity(step, *enerd, *inputrec);
2117 /* In case we don't have constraints and are using GPUs, the next balancing
2118 * region starts here.
2119 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
2120 * virial calculation and COM pulling, is not thus not included in
2121 * the balance timing, which is ok as most tasks do communication.
2123 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);