2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013-2019,2020,2021, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
49 #include "gromacs/applied_forces/awh/awh.h"
50 #include "gromacs/domdec/dlbtiming.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/domdec/gpuhaloexchange.h"
54 #include "gromacs/domdec/partition.h"
55 #include "gromacs/essentialdynamics/edsam.h"
56 #include "gromacs/ewald/pme.h"
57 #include "gromacs/ewald/pme_pp.h"
58 #include "gromacs/ewald/pme_pp_comm_gpu.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
61 #include "gromacs/gmxlib/nonbonded/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] += pull_potential(
207 pull_work, mdatoms->massT, &pbc, cr, t, lambda[efptRESTRAINT], as_rvec_array(x.data()), force, &dvdl);
208 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
209 wallcycle_stop(wcycle, ewcPULLPOT);
212 static void pme_receive_force_ener(t_forcerec* fr,
214 gmx::ForceWithVirial* forceWithVirial,
215 gmx_enerdata_t* enerd,
216 bool useGpuPmePpComms,
217 bool receivePmeForceToGpu,
218 gmx_wallcycle_t wcycle)
220 real e_q, e_lj, dvdl_q, dvdl_lj;
221 float cycles_ppdpme, cycles_seppme;
223 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
224 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
226 /* In case of node-splitting, the PP nodes receive the long-range
227 * forces, virial and energy from the PME nodes here.
229 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
232 gmx_pme_receive_f(fr->pmePpCommGpu.get(),
240 receivePmeForceToGpu,
242 enerd->term[F_COUL_RECIP] += e_q;
243 enerd->term[F_LJ_RECIP] += e_lj;
244 enerd->dvdl_lin[efptCOUL] += dvdl_q;
245 enerd->dvdl_lin[efptVDW] += dvdl_lj;
249 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
251 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
254 static void print_large_forces(FILE* fp,
259 ArrayRef<const RVec> x,
260 ArrayRef<const RVec> f)
262 real force2Tolerance = gmx::square(forceTolerance);
263 gmx::index numNonFinite = 0;
264 for (int i = 0; i < md->homenr; i++)
266 real force2 = norm2(f[i]);
267 bool nonFinite = !std::isfinite(force2);
268 if (force2 >= force2Tolerance || nonFinite)
271 "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
284 if (numNonFinite > 0)
286 /* Note that with MPI this fatal call on one rank might interrupt
287 * the printing on other ranks. But we can only avoid that with
288 * an expensive MPI barrier that we would need at each step.
290 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
294 //! When necessary, spreads forces on vsites and computes the virial for \p forceOutputs->forceWithShiftForces()
295 static void postProcessForceWithShiftForces(t_nrnb* nrnb,
296 gmx_wallcycle_t wcycle,
298 ArrayRef<const RVec> x,
299 ForceOutputs* forceOutputs,
301 const t_mdatoms& mdatoms,
302 const t_forcerec& fr,
303 gmx::VirtualSitesHandler* vsite,
304 const StepWorkload& stepWork)
306 ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
308 /* If we have NoVirSum forces, but we do not calculate the virial,
309 * we later sum the forceWithShiftForces buffer together with
310 * the noVirSum buffer and spread the combined vsite forces at once.
312 if (vsite && (!forceOutputs->haveForceWithVirial() || stepWork.computeVirial))
314 using VirialHandling = gmx::VirtualSitesHandler::VirialHandling;
316 auto f = forceWithShiftForces.force();
317 auto fshift = forceWithShiftForces.shiftForces();
318 const VirialHandling virialHandling =
319 (stepWork.computeVirial ? VirialHandling::Pbc : VirialHandling::None);
320 vsite->spreadForces(x, f, virialHandling, fshift, nullptr, nrnb, box, wcycle);
321 forceWithShiftForces.haveSpreadVsiteForces() = true;
324 if (stepWork.computeVirial)
326 /* Calculation of the virial must be done after vsites! */
328 0, mdatoms.homenr, as_rvec_array(x.data()), forceWithShiftForces, vir_force, box, nrnb, &fr, fr.pbcType);
332 //! Spread, compute virial for and sum forces, when necessary
333 static void postProcessForces(const t_commrec* cr,
336 gmx_wallcycle_t wcycle,
338 ArrayRef<const RVec> x,
339 ForceOutputs* forceOutputs,
341 const t_mdatoms* mdatoms,
342 const t_forcerec* fr,
343 gmx::VirtualSitesHandler* vsite,
344 const StepWorkload& stepWork)
346 // Extract the final output force buffer, which is also the buffer for forces with shift forces
347 ArrayRef<RVec> f = forceOutputs->forceWithShiftForces().force();
349 if (forceOutputs->haveForceWithVirial())
351 auto& forceWithVirial = forceOutputs->forceWithVirial();
355 /* Spread the mesh force on virtual sites to the other particles...
356 * This is parallellized. MPI communication is performed
357 * if the constructing atoms aren't local.
359 GMX_ASSERT(!stepWork.computeVirial || f.data() != forceWithVirial.force_.data(),
360 "We need separate force buffers for shift and virial forces when "
361 "computing the virial");
362 GMX_ASSERT(!stepWork.computeVirial
363 || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
364 "We should spread the force with shift forces separately when computing "
366 const gmx::VirtualSitesHandler::VirialHandling virialHandling =
367 (stepWork.computeVirial ? gmx::VirtualSitesHandler::VirialHandling::NonLinear
368 : gmx::VirtualSitesHandler::VirialHandling::None);
369 matrix virial = { { 0 } };
370 vsite->spreadForces(x, forceWithVirial.force_, virialHandling, {}, virial, nrnb, box, wcycle);
371 forceWithVirial.addVirialContribution(virial);
374 if (stepWork.computeVirial)
376 /* Now add the forces, this is local */
377 sum_forces(f, forceWithVirial.force_);
379 /* Add the direct virial contributions */
381 forceWithVirial.computeVirial_,
382 "forceWithVirial should request virial computation when we request the virial");
383 m_add(vir_force, forceWithVirial.getVirial(), vir_force);
387 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
393 GMX_ASSERT(vsite == nullptr || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
394 "We should have spread the vsite forces (earlier)");
397 if (fr->print_force >= 0)
399 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
403 static void do_nb_verlet(t_forcerec* fr,
404 const interaction_const_t* ic,
405 gmx_enerdata_t* enerd,
406 const StepWorkload& stepWork,
407 const InteractionLocality ilocality,
411 gmx_wallcycle_t wcycle)
413 if (!stepWork.computeNonbondedForces)
415 /* skip non-bonded calculation */
419 nonbonded_verlet_t* nbv = fr->nbv.get();
421 /* GPU kernel launch overhead is already timed separately */
424 /* When dynamic pair-list pruning is requested, we need to prune
425 * at nstlistPrune steps.
427 if (nbv->isDynamicPruningStepCpu(step))
429 /* Prune the pair-list beyond fr->ic->rlistPrune using
430 * the current coordinates of the atoms.
432 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
433 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
434 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
438 nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
441 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
443 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, v.ssize());
445 /* Note that we would like to avoid this conditional by putting it
446 * into the omp pragma instead, but then we still take the full
447 * omp parallel for overhead (at least with gcc5).
449 if (!useOpenmpThreading || nth == 1)
458 #pragma omp parallel for num_threads(nth) schedule(static)
459 for (gmx::index i = 0; i < v.ssize(); i++)
466 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
468 * \param groupOptions Group options, containing T-coupling options
470 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
472 real nrdfCoupled = 0;
473 real nrdfUncoupled = 0;
474 real kineticEnergy = 0;
475 for (int g = 0; g < groupOptions.ngtc; g++)
477 if (groupOptions.tau_t[g] >= 0)
479 nrdfCoupled += groupOptions.nrdf[g];
480 kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * BOLTZ;
484 nrdfUncoupled += groupOptions.nrdf[g];
488 /* This conditional with > also catches nrdf=0 */
489 if (nrdfCoupled > nrdfUncoupled)
491 return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
499 /*! \brief This routine checks that the potential energy is finite.
501 * Always checks that the potential energy is finite. If step equals
502 * inputrec.init_step also checks that the magnitude of the potential energy
503 * is reasonable. Terminates with a fatal error when a check fails.
504 * Note that passing this check does not guarantee finite forces,
505 * since those use slightly different arithmetics. But in most cases
506 * there is just a narrow coordinate range where forces are not finite
507 * and energies are finite.
509 * \param[in] step The step number, used for checking and printing
510 * \param[in] enerd The energy data; the non-bonded group energies need to be added to
511 * enerd.term[F_EPOT] before calling this routine \param[in] inputrec The input record
513 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
515 /* Threshold valid for comparing absolute potential energy against
516 * the kinetic energy. Normally one should not consider absolute
517 * potential energy values, but with a factor of one million
518 * we should never get false positives.
520 constexpr real c_thresholdFactor = 1e6;
522 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
523 real averageKineticEnergy = 0;
524 /* We only check for large potential energy at the initial step,
525 * because that is by far the most likely step for this too occur
526 * and because computing the average kinetic energy is not free.
527 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
528 * before they become NaN.
530 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
532 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
535 if (energyIsNotFinite
536 || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
541 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
542 "contributions to the energy are %g and %g, respectively. A %s potential energy "
543 "can be caused by overlapping interactions in bonded interactions or very large%s "
544 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
545 "configuration, incorrect interactions or parameters in the topology.",
548 energyIsNotFinite ? "not finite" : "extremely high",
550 enerd.term[F_COUL_SR],
551 energyIsNotFinite ? "non-finite" : "very high",
552 energyIsNotFinite ? " or Nan" : "");
556 /*! \brief Return true if there are special forces computed this step.
558 * The conditionals exactly correspond to those in computeSpecialForces().
560 static bool haveSpecialForces(const t_inputrec& inputrec,
561 const gmx::ForceProviders& forceProviders,
562 const pull_t* pull_work,
563 const bool computeForces,
567 return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
568 (inputrec.bPull && pull_have_potential(*pull_work)) || // pull
569 inputrec.bRot || // enforced rotation
570 (ed != nullptr) || // flooding
571 (inputrec.bIMD && computeForces)); // IMD
574 /*! \brief Compute forces and/or energies for special algorithms
576 * The intention is to collect all calls to algorithms that compute
577 * forces on local atoms only and that do not contribute to the local
578 * virial sum (but add their virial contribution separately).
579 * Eventually these should likely all become ForceProviders.
580 * Within this function the intention is to have algorithms that do
581 * global communication at the end, so global barriers within the MD loop
582 * are as close together as possible.
584 * \param[in] fplog The log file
585 * \param[in] cr The communication record
586 * \param[in] inputrec The input record
587 * \param[in] awh The Awh module (nullptr if none in use).
588 * \param[in] enforcedRotation Enforced rotation module.
589 * \param[in] imdSession The IMD session
590 * \param[in] pull_work The pull work structure.
591 * \param[in] step The current MD step
592 * \param[in] t The current time
593 * \param[in,out] wcycle Wallcycle accounting struct
594 * \param[in,out] forceProviders Pointer to a list of force providers
595 * \param[in] box The unit cell
596 * \param[in] x The coordinates
597 * \param[in] mdatoms Per atom properties
598 * \param[in] lambda Array of free-energy lambda values
599 * \param[in] stepWork Step schedule flags
600 * \param[in,out] forceWithVirialMtsLevel0 Force and virial for MTS level0 forces
601 * \param[in,out] forceWithVirialMtsLevel1 Force and virial for MTS level1 forces, can be nullptr
602 * \param[in,out] enerd Energy buffer
603 * \param[in,out] ed Essential dynamics pointer
604 * \param[in] didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
606 * \todo Remove didNeighborSearch, which is used incorrectly.
607 * \todo Convert all other algorithms called here to ForceProviders.
609 static void computeSpecialForces(FILE* fplog,
611 const t_inputrec* inputrec,
613 gmx_enfrot* enforcedRotation,
614 gmx::ImdSession* imdSession,
618 gmx_wallcycle_t wcycle,
619 gmx::ForceProviders* forceProviders,
621 gmx::ArrayRef<const gmx::RVec> x,
622 const t_mdatoms* mdatoms,
623 gmx::ArrayRef<const real> lambda,
624 const StepWorkload& stepWork,
625 gmx::ForceWithVirial* forceWithVirialMtsLevel0,
626 gmx::ForceWithVirial* forceWithVirialMtsLevel1,
627 gmx_enerdata_t* enerd,
629 bool didNeighborSearch)
631 /* NOTE: Currently all ForceProviders only provide forces.
632 * When they also provide energies, remove this conditional.
634 if (stepWork.computeForces)
636 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
637 gmx::ForceProviderOutput forceProviderOutput(forceWithVirialMtsLevel0, enerd);
639 /* Collect forces from modules */
640 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
643 if (inputrec->bPull && pull_have_potential(*pull_work))
645 const int mtsLevel = forceGroupMtsLevel(inputrec->mtsLevels, gmx::MtsForceGroups::Pull);
646 if (mtsLevel == 0 || stepWork.computeSlowForces)
648 auto& forceWithVirial = (mtsLevel == 0) ? forceWithVirialMtsLevel0 : forceWithVirialMtsLevel1;
649 pull_potential_wrapper(
650 cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work, lambda.data(), t, wcycle);
655 const int mtsLevel = forceGroupMtsLevel(inputrec->mtsLevels, gmx::MtsForceGroups::Pull);
656 if (mtsLevel == 0 || stepWork.computeSlowForces)
658 const bool needForeignEnergyDifferences = awh->needForeignEnergyDifferences(step);
659 std::vector<double> foreignLambdaDeltaH, foreignLambdaDhDl;
660 if (needForeignEnergyDifferences)
662 enerd->foreignLambdaTerms.finalizePotentialContributions(
663 enerd->dvdl_lin, lambda, *inputrec->fepvals);
664 std::tie(foreignLambdaDeltaH, foreignLambdaDhDl) = enerd->foreignLambdaTerms.getTerms(cr);
667 auto& forceWithVirial = (mtsLevel == 0) ? forceWithVirialMtsLevel0 : forceWithVirialMtsLevel1;
668 enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(inputrec->pbcType,
681 rvec* f = as_rvec_array(forceWithVirialMtsLevel0->force_.data());
683 /* Add the forces from enforced rotation potentials (if any) */
686 wallcycle_start(wcycle, ewcROTadd);
687 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
688 wallcycle_stop(wcycle, ewcROTadd);
693 /* Note that since init_edsam() is called after the initialization
694 * of forcerec, edsam doesn't request the noVirSum force buffer.
695 * Thus if no other algorithm (e.g. PME) requires it, the forces
696 * here will contribute to the virial.
698 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
701 /* Add forces from interactive molecular dynamics (IMD), if any */
702 if (inputrec->bIMD && stepWork.computeForces)
704 imdSession->applyForces(f);
708 /*! \brief Launch the prepare_step and spread stages of PME GPU.
710 * \param[in] pmedata The PME structure
711 * \param[in] box The box matrix
712 * \param[in] stepWork Step schedule flags
713 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in the device memory.
714 * \param[in] lambdaQ The Coulomb lambda of the current state.
715 * \param[in] wcycle The wallcycle structure
717 static inline void launchPmeGpuSpread(gmx_pme_t* pmedata,
719 const StepWorkload& stepWork,
720 GpuEventSynchronizer* xReadyOnDevice,
722 gmx_wallcycle_t wcycle)
724 pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
725 pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle, lambdaQ);
728 /*! \brief Launch the FFT and gather stages of PME GPU
730 * This function only implements setting the output forces (no accumulation).
732 * \param[in] pmedata The PME structure
733 * \param[in] lambdaQ The Coulomb lambda of the current system state.
734 * \param[in] wcycle The wallcycle structure
735 * \param[in] stepWork Step schedule flags
737 static void launchPmeGpuFftAndGather(gmx_pme_t* pmedata,
739 gmx_wallcycle_t wcycle,
740 const gmx::StepWorkload& stepWork)
742 pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
743 pme_gpu_launch_gather(pmedata, wcycle, lambdaQ);
747 * Polling wait for either of the PME or nonbonded GPU tasks.
749 * Instead of a static order in waiting for GPU tasks, this function
750 * polls checking which of the two tasks completes first, and does the
751 * associated force buffer reduction overlapped with the other task.
752 * By doing that, unlike static scheduling order, it can always overlap
753 * one of the reductions, regardless of the GPU task completion order.
755 * \param[in] nbv Nonbonded verlet structure
756 * \param[in,out] pmedata PME module data
757 * \param[in,out] forceOutputsNonbonded Force outputs for the non-bonded forces and shift forces
758 * \param[in,out] forceOutputsPme Force outputs for the PME forces and virial
759 * \param[in,out] enerd Energy data structure results are reduced into
760 * \param[in] lambdaQ The Coulomb lambda of the current system state.
761 * \param[in] stepWork Step schedule flags
762 * \param[in] wcycle The wallcycle structure
764 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
766 gmx::ForceOutputs* forceOutputsNonbonded,
767 gmx::ForceOutputs* forceOutputsPme,
768 gmx_enerdata_t* enerd,
770 const StepWorkload& stepWork,
771 gmx_wallcycle_t wcycle)
773 bool isPmeGpuDone = false;
774 bool isNbGpuDone = false;
776 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
778 while (!isPmeGpuDone || !isNbGpuDone)
782 GpuTaskCompletion completionType =
783 (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
784 isPmeGpuDone = pme_gpu_try_finish_task(
785 pmedata, stepWork, wcycle, &forceOutputsPme->forceWithVirial(), enerd, lambdaQ, completionType);
790 auto& forceBuffersNonbonded = forceOutputsNonbonded->forceWithShiftForces();
791 GpuTaskCompletion completionType =
792 (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
793 isNbGpuDone = Nbnxm::gpu_try_finish_task(nbv->gpu_nbv,
796 enerd->grpp.ener[egLJSR].data(),
797 enerd->grpp.ener[egCOULSR].data(),
798 forceBuffersNonbonded.shiftForces(),
804 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceBuffersNonbonded.force());
810 /*! \brief Set up the different force buffers; also does clearing.
812 * \param[in] forceHelperBuffers Helper force buffers
813 * \param[in] force force array
814 * \param[in] stepWork Step schedule flags
815 * \param[out] wcycle wallcycle recording structure
817 * \returns Cleared force output structure
819 static ForceOutputs setupForceOutputs(ForceHelperBuffers* forceHelperBuffers,
820 gmx::ArrayRefWithPadding<gmx::RVec> force,
821 const StepWorkload& stepWork,
822 gmx_wallcycle_t wcycle)
824 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
826 /* NOTE: We assume fr->shiftForces is all zeros here */
827 gmx::ForceWithShiftForces forceWithShiftForces(
828 force, stepWork.computeVirial, forceHelperBuffers->shiftForces());
830 if (stepWork.computeForces)
832 /* Clear the short- and long-range forces */
833 clearRVecs(forceWithShiftForces.force(), true);
835 /* Clear the shift forces */
836 clearRVecs(forceWithShiftForces.shiftForces(), false);
839 /* If we need to compute the virial, we might need a separate
840 * force buffer for algorithms for which the virial is calculated
841 * directly, such as PME. Otherwise, forceWithVirial uses the
842 * the same force (f in legacy calls) buffer as other algorithms.
844 const bool useSeparateForceWithVirialBuffer =
845 (stepWork.computeForces
846 && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
847 /* forceWithVirial uses the local atom range only */
848 gmx::ForceWithVirial forceWithVirial(
849 useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
850 : force.unpaddedArrayRef(),
851 stepWork.computeVirial);
853 if (useSeparateForceWithVirialBuffer)
855 /* TODO: update comment
856 * We only compute forces on local atoms. Note that vsites can
857 * spread to non-local atoms, but that part of the buffer is
858 * cleared separately in the vsite spreading code.
860 clearRVecs(forceWithVirial.force_, true);
863 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
866 forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(), forceWithVirial);
870 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
872 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec& inputrec,
873 const t_forcerec& fr,
874 const pull_t* pull_work,
876 const t_mdatoms& mdatoms,
877 const SimulationWorkload& simulationWork,
878 const StepWorkload& stepWork)
880 DomainLifetimeWorkload domainWork;
881 // Note that haveSpecialForces is constant over the whole run
882 domainWork.haveSpecialForces =
883 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
884 domainWork.haveCpuListedForceWork = false;
885 domainWork.haveCpuBondedWork = false;
886 for (const auto& listedForces : fr.listedForces)
888 if (listedForces.haveCpuListedForces(*fr.fcdata))
890 domainWork.haveCpuListedForceWork = true;
892 if (listedForces.haveCpuBondeds())
894 domainWork.haveCpuBondedWork = true;
897 domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
898 // Note that haveFreeEnergyWork is constant over the whole run
899 domainWork.haveFreeEnergyWork = (fr.efep != efepNO && mdatoms.nPerturbed != 0);
900 // We assume we have local force work if there are CPU
901 // force tasks including PME or nonbondeds.
902 domainWork.haveCpuLocalForceWork =
903 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
904 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
905 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
910 /*! \brief Set up force flag stuct from the force bitmask.
912 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
913 * \param[in] mtsLevels The multiple time-stepping levels, either empty or 2 levels
914 * \param[in] step The current MD step
915 * \param[in] simulationWork Simulation workload description.
916 * \param[in] rankHasPmeDuty If this rank computes PME.
918 * \returns New Stepworkload description.
920 static StepWorkload setupStepWorkload(const int legacyFlags,
921 ArrayRef<const gmx::MtsLevel> mtsLevels,
923 const SimulationWorkload& simulationWork,
924 const bool rankHasPmeDuty)
926 GMX_ASSERT(mtsLevels.empty() || mtsLevels.size() == 2, "Expect 0 or 2 MTS levels");
927 const bool computeSlowForces = (mtsLevels.empty() || step % mtsLevels[1].stepFactor == 0);
930 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
931 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
932 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
933 flags.computeSlowForces = computeSlowForces;
934 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
935 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
936 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
937 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
938 flags.computeNonbondedForces =
939 ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded
940 && !(simulationWork.computeNonbondedAtMtsLevel1 && !computeSlowForces);
941 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
943 if (simulationWork.useGpuBufferOps)
945 GMX_ASSERT(simulationWork.useGpuNonbonded,
946 "Can only offload buffer ops if nonbonded computation is also offloaded");
948 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
949 // on virial steps the CPU reduction path is taken
950 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
951 flags.useGpuPmeFReduction = flags.computeSlowForces && flags.useGpuFBufferOps && simulationWork.useGpuPme
952 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication);
953 flags.useGpuXHalo = simulationWork.useGpuHaloExchange;
954 flags.useGpuFHalo = simulationWork.useGpuHaloExchange && flags.useGpuFBufferOps;
960 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
962 * TODO: eliminate \p useGpuPmeOnThisRank when this is
963 * incorporated in DomainLifetimeWorkload.
965 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
966 gmx::GpuBonded* gpuBonded,
968 gmx_enerdata_t* enerd,
969 const gmx::MdrunScheduleWorkload& runScheduleWork,
970 bool useGpuPmeOnThisRank,
972 gmx_wallcycle_t wcycle)
974 if (runScheduleWork.simulationWork.useGpuNonbonded && runScheduleWork.stepWork.computeNonbondedForces)
976 /* Launch pruning before buffer clearing because the API overhead of the
977 * clear kernel launches can leave the GPU idle while it could be running
980 if (nbv->isDynamicPruningStepGpu(step))
982 nbv->dispatchPruneKernelGpu(step);
985 /* now clear the GPU outputs while we finish the step on the CPU */
986 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
987 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
988 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
989 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
990 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
993 if (useGpuPmeOnThisRank)
995 pme_gpu_reinit_computation(pmedata, wcycle);
998 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
1000 // in principle this should be included in the DD balancing region,
1001 // but generally it is infrequent so we'll omit it for the sake of
1003 gpuBonded->waitAccumulateEnergyTerms(enerd);
1005 gpuBonded->clearEnergies();
1009 //! \brief Data structure to hold dipole-related data and staging arrays
1012 //! Dipole staging for fast summing over MPI
1013 gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
1014 //! Dipole staging for states A and B (index 0 and 1 resp.)
1015 gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
1019 static void reduceAndUpdateMuTot(DipoleData* dipoleData,
1020 const t_commrec* cr,
1021 const bool haveFreeEnergy,
1022 gmx::ArrayRef<const real> lambda,
1024 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1028 gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
1029 ddBalanceRegionHandler.reopenRegionCpu();
1031 for (int i = 0; i < 2; i++)
1033 for (int j = 0; j < DIM; j++)
1035 dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
1039 if (!haveFreeEnergy)
1041 copy_rvec(dipoleData->muStateAB[0], muTotal);
1045 for (int j = 0; j < DIM; j++)
1047 muTotal[j] = (1.0 - lambda[efptCOUL]) * dipoleData->muStateAB[0][j]
1048 + lambda[efptCOUL] * dipoleData->muStateAB[1][j];
1053 /*! \brief Combines MTS level0 and level1 force buffes into a full and MTS-combined force buffer.
1055 * \param[in] numAtoms The number of atoms to combine forces for
1056 * \param[in,out] forceMtsLevel0 Input: F_level0, output: F_level0 + F_level1
1057 * \param[in,out] forceMts Input: F_level1, output: F_level0 + mtsFactor * F_level1
1058 * \param[in] mtsFactor The factor between the level0 and level1 time step
1060 static void combineMtsForces(const int numAtoms,
1061 ArrayRef<RVec> forceMtsLevel0,
1062 ArrayRef<RVec> forceMts,
1063 const real mtsFactor)
1065 const int gmx_unused numThreads = gmx_omp_nthreads_get(emntDefault);
1066 #pragma omp parallel for num_threads(numThreads) schedule(static)
1067 for (int i = 0; i < numAtoms; i++)
1069 const RVec forceMtsLevel0Tmp = forceMtsLevel0[i];
1070 forceMtsLevel0[i] += forceMts[i];
1071 forceMts[i] = forceMtsLevel0Tmp + mtsFactor * forceMts[i];
1075 /*! \brief Setup for the local and non-local GPU force reductions:
1076 * reinitialization plus the registration of forces and dependencies.
1078 * \param [in] runScheduleWork Schedule workload flag structure
1079 * \param [in] cr Communication record object
1080 * \param [in] fr Force record object
1082 static void setupGpuForceReductions(gmx::MdrunScheduleWorkload* runScheduleWork,
1083 const t_commrec* cr,
1087 nonbonded_verlet_t* nbv = fr->nbv.get();
1088 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1090 // (re-)initialize local GPU force reduction
1091 const bool accumulate =
1092 runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr);
1093 const int atomStart = 0;
1094 fr->gpuForceReduction[gmx::AtomLocality::Local]->reinit(stateGpu->getForces(),
1095 nbv->getNumAtoms(AtomLocality::Local),
1096 nbv->getGridIndices(),
1099 stateGpu->fReducedOnDevice());
1101 // register forces and add dependencies
1102 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerNbnxmForce(nbv->getGpuForces());
1104 if (runScheduleWork->simulationWork.useGpuPme
1105 && (thisRankHasDuty(cr, DUTY_PME) || runScheduleWork->simulationWork.useGpuPmePpCommunication))
1107 void* forcePtr = thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1108 : // PME force buffer on same GPU
1109 fr->pmePpCommGpu->getGpuForceStagingPtr(); // buffer received from other GPU
1110 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerRvecForce(forcePtr);
1112 GpuEventSynchronizer* const pmeSynchronizer =
1113 (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1114 : // PME force buffer on same GPU
1115 fr->pmePpCommGpu->getForcesReadySynchronizer()); // buffer received from other GPU
1116 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(pmeSynchronizer);
1119 if ((runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr))
1120 && !runScheduleWork->simulationWork.useGpuHaloExchange)
1122 auto forcesReadyLocality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
1123 const bool useGpuForceBufferOps = true;
1124 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1125 stateGpu->getForcesReadyOnDeviceEvent(forcesReadyLocality, useGpuForceBufferOps));
1128 if (runScheduleWork->simulationWork.useGpuHaloExchange)
1130 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1131 cr->dd->gpuHaloExchange[0][0]->getForcesReadyOnDeviceEvent());
1134 if (havePPDomainDecomposition(cr))
1136 // (re-)initialize non-local GPU force reduction
1137 const bool accumulate = runScheduleWork->domainWork.haveCpuBondedWork
1138 || runScheduleWork->domainWork.haveFreeEnergyWork;
1139 const int atomStart = dd_numHomeAtoms(*cr->dd);
1140 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->reinit(stateGpu->getForces(),
1141 nbv->getNumAtoms(AtomLocality::NonLocal),
1142 nbv->getGridIndices(),
1146 // register forces and add dependencies
1147 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->registerNbnxmForce(nbv->getGpuForces());
1148 if (runScheduleWork->domainWork.haveCpuBondedWork || runScheduleWork->domainWork.haveFreeEnergyWork)
1150 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->addDependency(
1151 stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::NonLocal, true));
1157 void do_force(FILE* fplog,
1158 const t_commrec* cr,
1159 const gmx_multisim_t* ms,
1160 const t_inputrec* inputrec,
1162 gmx_enfrot* enforcedRotation,
1163 gmx::ImdSession* imdSession,
1167 gmx_wallcycle_t wcycle,
1168 const gmx_localtop_t* top,
1170 gmx::ArrayRefWithPadding<gmx::RVec> x,
1171 const history_t* hist,
1172 gmx::ForceBuffersView* forceView,
1174 const t_mdatoms* mdatoms,
1175 gmx_enerdata_t* enerd,
1176 gmx::ArrayRef<const real> lambda,
1178 gmx::MdrunScheduleWorkload* runScheduleWork,
1179 gmx::VirtualSitesHandler* vsite,
1184 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1186 auto force = forceView->forceWithPadding();
1187 GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1188 "The size of the force buffer should be at least the number of atoms to compute "
1191 nonbonded_verlet_t* nbv = fr->nbv.get();
1192 interaction_const_t* ic = fr->ic;
1194 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1196 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1198 runScheduleWork->stepWork = setupStepWorkload(
1199 legacyFlags, inputrec->mtsLevels, step, simulationWork, thisRankHasDuty(cr, DUTY_PME));
1200 const StepWorkload& stepWork = runScheduleWork->stepWork;
1202 const bool useGpuPmeOnThisRank =
1203 simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces;
1205 /* At a search step we need to start the first balancing region
1206 * somewhere early inside the step after communication during domain
1207 * decomposition (and not during the previous step as usual).
1209 if (stepWork.doNeighborSearch)
1211 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1214 clear_mat(vir_force);
1216 if (fr->pbcType != PbcType::No)
1218 /* Compute shift vectors every step,
1219 * because of pressure coupling or box deformation!
1221 if (stepWork.haveDynamicBox && stepWork.stateChanged)
1223 calc_shifts(box, fr->shift_vec);
1226 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1227 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1230 put_atoms_in_box_omp(fr->pbcType,
1232 x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1233 gmx_omp_nthreads_get(emntDefault));
1234 inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1238 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1240 const bool pmeSendCoordinatesFromGpu =
1241 GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1242 const bool reinitGpuPmePpComms =
1243 GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1245 const auto localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1246 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1247 AtomLocality::Local, simulationWork, stepWork)
1250 // Copy coordinate from the GPU if update is on the GPU and there
1251 // are forces to be computed on the CPU, or for the computation of
1252 // virial, or if host-side data will be transferred from this task
1253 // to a remote task for halo exchange or PME-PP communication. At
1254 // search steps the current coordinates are already on the host,
1255 // hence copy is not needed.
1256 const bool haveHostPmePpComms =
1257 !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1259 GMX_ASSERT(simulationWork.useGpuHaloExchange
1260 == ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange[0].empty())),
1261 "The GPU halo exchange is active, but it has not been constructed.");
1262 const bool haveHostHaloExchangeComms =
1263 havePPDomainDecomposition(cr) && !simulationWork.useGpuHaloExchange;
1265 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1266 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1267 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1268 || haveHostPmePpComms || haveHostHaloExchangeComms))
1270 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1271 haveCopiedXFromGpu = true;
1274 // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1275 // Otherwise the send will occur after H2D coordinate transfer.
1276 if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu && stepWork.computeSlowForces)
1278 /* Send particle coordinates to the pme nodes */
1279 if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1281 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1284 gmx_pme_send_coordinates(fr,
1287 as_rvec_array(x.unpaddedArrayRef().data()),
1290 (stepWork.computeVirial || stepWork.computeEnergy),
1292 simulationWork.useGpuPmePpCommunication,
1293 reinitGpuPmePpComms,
1294 pmeSendCoordinatesFromGpu,
1295 localXReadyOnDevice,
1299 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1300 // The local coordinates can be copied right away.
1301 // NOTE: Consider moving this copy to right after they are updated and constrained,
1302 // if the later is not offloaded.
1303 if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1305 if (stepWork.doNeighborSearch)
1307 // TODO refactor this to do_md, after partitioning.
1308 stateGpu->reinit(mdatoms->homenr,
1309 cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1310 if (useGpuPmeOnThisRank)
1312 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1313 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1316 // We need to copy coordinates when:
1317 // 1. Update is not offloaded
1318 // 2. The buffers were reinitialized on search step
1319 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1321 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1322 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1326 // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1327 // Otherwise the send will occur before the H2D coordinate transfer.
1328 if (!thisRankHasDuty(cr, DUTY_PME) && pmeSendCoordinatesFromGpu)
1330 /* Send particle coordinates to the pme nodes */
1331 gmx_pme_send_coordinates(fr,
1334 as_rvec_array(x.unpaddedArrayRef().data()),
1337 (stepWork.computeVirial || stepWork.computeEnergy),
1339 simulationWork.useGpuPmePpCommunication,
1340 reinitGpuPmePpComms,
1341 pmeSendCoordinatesFromGpu,
1342 localXReadyOnDevice,
1346 if (useGpuPmeOnThisRank)
1348 launchPmeGpuSpread(fr->pmedata, box, stepWork, localXReadyOnDevice, lambda[efptCOUL], wcycle);
1351 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1353 /* do gridding for pair search */
1354 if (stepWork.doNeighborSearch)
1356 if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1358 fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1362 // - vzero is constant, do we need to pass it?
1363 // - box_diag should be passed directly to nbnxn_put_on_grid
1369 box_diag[XX] = box[XX][XX];
1370 box_diag[YY] = box[YY][YY];
1371 box_diag[ZZ] = box[ZZ][ZZ];
1373 wallcycle_start(wcycle, ewcNS);
1374 if (!DOMAINDECOMP(cr))
1376 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1377 nbnxn_put_on_grid(nbv,
1383 { 0, mdatoms->homenr },
1386 x.unpaddedArrayRef(),
1389 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1393 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1394 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1395 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1398 nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1399 gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1402 wallcycle_stop(wcycle, ewcNS);
1404 /* initialize the GPU nbnxm atom data and bonded data structures */
1405 if (simulationWork.useGpuNonbonded)
1407 // Note: cycle counting only nononbondeds, gpuBonded counts internally
1408 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1409 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1410 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1411 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1412 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1416 /* Now we put all atoms on the grid, we can assign bonded
1417 * interactions to the GPU, where the grid order is
1418 * needed. Also the xq, f and fshift device buffers have
1419 * been reallocated if needed, so the bonded code can
1420 * learn about them. */
1421 // TODO the xq, f, and fshift buffers are now shared
1422 // resources, so they should be maintained by a
1423 // higher-level object than the nb module.
1424 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1426 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1427 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1428 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1432 // Need to run after the GPU-offload bonded interaction lists
1433 // are set up to be able to determine whether there is bonded work.
1434 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1435 *inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1437 wallcycle_start_nocount(wcycle, ewcNS);
1438 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1439 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1440 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1442 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1444 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1445 wallcycle_stop(wcycle, ewcNS);
1447 if (stepWork.useGpuXBufferOps)
1449 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1452 if (simulationWork.useGpuBufferOps)
1454 setupGpuForceReductions(runScheduleWork, cr, fr);
1457 else if (!EI_TPI(inputrec->eI) && stepWork.computeNonbondedForces)
1459 if (stepWork.useGpuXBufferOps)
1461 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1462 nbv->convertCoordinatesGpu(
1463 AtomLocality::Local, false, stateGpu->getCoordinates(), localXReadyOnDevice);
1467 if (simulationWork.useGpuUpdate)
1469 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1470 GMX_ASSERT(haveCopiedXFromGpu,
1471 "a wait should only be triggered if copy has been scheduled");
1472 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1474 nbv->convertCoordinates(AtomLocality::Local, false, x.unpaddedArrayRef());
1478 if (simulationWork.useGpuNonbonded && (stepWork.computeNonbondedForces || domainWork.haveGpuBondedWork))
1480 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1482 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1483 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1484 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1485 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1487 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1489 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1490 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1491 // with X buffer ops offloaded to the GPU on all but the search steps
1493 // bonded work not split into separate local and non-local, so with DD
1494 // we can only launch the kernel after non-local coordinates have been received.
1495 if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1497 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1500 /* launch local nonbonded work on GPU */
1501 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1502 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1503 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1504 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1505 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1508 if (useGpuPmeOnThisRank)
1510 // In PME GPU and mixed mode we launch FFT / gather after the
1511 // X copy/transform to allow overlap as well as after the GPU NB
1512 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1513 // the nonbonded kernel.
1514 launchPmeGpuFftAndGather(fr->pmedata, lambda[efptCOUL], wcycle, stepWork);
1517 /* Communicate coordinates and sum dipole if necessary +
1518 do non-local pair search */
1519 if (havePPDomainDecomposition(cr))
1521 if (stepWork.doNeighborSearch)
1523 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1524 wallcycle_start_nocount(wcycle, ewcNS);
1525 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1526 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1527 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1529 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1530 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1531 wallcycle_stop(wcycle, ewcNS);
1532 // TODO refactor this GPU halo exchange re-initialisation
1533 // to location in do_md where GPU halo exchange is
1534 // constructed at partitioning, after above stateGpu
1535 // re-initialization has similarly been refactored
1536 if (simulationWork.useGpuHaloExchange)
1538 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1543 if (stepWork.useGpuXHalo)
1545 // The following must be called after local setCoordinates (which records an event
1546 // when the coordinate data has been copied to the device).
1547 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1549 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1551 // non-local part of coordinate buffer must be copied back to host for CPU work
1552 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1557 if (simulationWork.useGpuUpdate)
1559 GMX_ASSERT(haveCopiedXFromGpu,
1560 "a wait should only be triggered if copy has been scheduled");
1561 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1563 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1566 if (stepWork.useGpuXBufferOps)
1568 if (!useGpuPmeOnThisRank && !stepWork.useGpuXHalo)
1570 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1572 nbv->convertCoordinatesGpu(AtomLocality::NonLocal,
1574 stateGpu->getCoordinates(),
1575 stateGpu->getCoordinatesReadyOnDeviceEvent(
1576 AtomLocality::NonLocal, simulationWork, stepWork));
1580 nbv->convertCoordinates(AtomLocality::NonLocal, false, x.unpaddedArrayRef());
1584 if (simulationWork.useGpuNonbonded)
1587 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1589 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1590 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1591 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1592 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1593 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1596 if (domainWork.haveGpuBondedWork)
1598 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1601 /* launch non-local nonbonded tasks on GPU */
1602 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1603 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1604 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1605 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1606 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1610 if (simulationWork.useGpuNonbonded && stepWork.computeNonbondedForces)
1612 /* launch D2H copy-back F */
1613 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1614 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1616 if (havePPDomainDecomposition(cr))
1618 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1620 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1621 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1623 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1625 fr->gpuBonded->launchEnergyTransfer();
1627 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1630 gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1631 if (fr->wholeMoleculeTransform)
1633 xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1636 DipoleData dipoleData;
1638 if (simulationWork.computeMuTot)
1640 const int start = 0;
1642 /* Calculate total (local) dipole moment in a temporary common array.
1643 * This makes it possible to sum them over nodes faster.
1645 gmx::ArrayRef<const gmx::RVec> xRef =
1646 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1652 mdatoms->nChargePerturbed,
1653 dipoleData.muStaging[0],
1654 dipoleData.muStaging[1]);
1656 reduceAndUpdateMuTot(&dipoleData, cr, (fr->efep != efepNO), lambda, muTotal, ddBalanceRegionHandler);
1659 /* Reset energies */
1660 reset_enerdata(enerd);
1662 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1664 wallcycle_start(wcycle, ewcPPDURINGPME);
1665 dd_force_flop_start(cr->dd, nrnb);
1668 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1669 // this wait ensures that the D2H transfer is complete.
1670 if ((simulationWork.useGpuUpdate)
1671 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1673 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1678 wallcycle_start(wcycle, ewcROT);
1679 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, stepWork.doNeighborSearch);
1680 wallcycle_stop(wcycle, ewcROT);
1683 /* Start the force cycle counter.
1684 * Note that a different counter is used for dynamic load balancing.
1686 wallcycle_start(wcycle, ewcFORCE);
1688 /* Set up and clear force outputs:
1689 * forceOutMtsLevel0: everything except what is in the other two outputs
1690 * forceOutMtsLevel1: PME-mesh and listed-forces group 1
1691 * forceOutNonbonded: non-bonded forces
1692 * Without multiple time stepping all point to the same object.
1693 * With multiple time-stepping the use is different for MTS fast (level0 only) and slow steps.
1695 ForceOutputs forceOutMtsLevel0 =
1696 setupForceOutputs(&fr->forceHelperBuffers[0], force, stepWork, wcycle);
1698 // Force output for MTS combined forces, only set at level1 MTS steps
1699 std::optional<ForceOutputs> forceOutMts =
1700 (fr->useMts && stepWork.computeSlowForces)
1701 ? std::optional(setupForceOutputs(&fr->forceHelperBuffers[1],
1702 forceView->forceMtsCombinedWithPadding(),
1707 ForceOutputs* forceOutMtsLevel1 =
1708 fr->useMts ? (stepWork.computeSlowForces ? &forceOutMts.value() : nullptr) : &forceOutMtsLevel0;
1710 const bool nonbondedAtMtsLevel1 = runScheduleWork->simulationWork.computeNonbondedAtMtsLevel1;
1712 ForceOutputs* forceOutNonbonded = nonbondedAtMtsLevel1 ? forceOutMtsLevel1 : &forceOutMtsLevel0;
1714 if (inputrec->bPull && pull_have_constraint(*pull_work))
1716 clear_pull_forces(pull_work);
1719 /* We calculate the non-bonded forces, when done on the CPU, here.
1720 * We do this before calling do_force_lowlevel, because in that
1721 * function, the listed forces are calculated before PME, which
1722 * does communication. With this order, non-bonded and listed
1723 * force calculation imbalance can be balanced out by the domain
1724 * decomposition load balancing.
1727 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1729 if (!useOrEmulateGpuNb)
1731 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1734 if (fr->efep != efepNO && stepWork.computeNonbondedForces)
1736 /* Calculate the local and non-local free energy interactions here.
1737 * Happens here on the CPU both with and without GPU.
1739 nbv->dispatchFreeEnergyKernel(InteractionLocality::Local,
1741 as_rvec_array(x.unpaddedArrayRef().data()),
1742 &forceOutNonbonded->forceWithShiftForces(),
1750 if (havePPDomainDecomposition(cr))
1752 nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal,
1754 as_rvec_array(x.unpaddedArrayRef().data()),
1755 &forceOutNonbonded->forceWithShiftForces(),
1765 if (stepWork.computeNonbondedForces && !useOrEmulateGpuNb)
1767 if (havePPDomainDecomposition(cr))
1769 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1772 if (stepWork.computeForces)
1774 /* Add all the non-bonded force to the normal force array.
1775 * This can be split into a local and a non-local part when overlapping
1776 * communication with calculation with domain decomposition.
1778 wallcycle_stop(wcycle, ewcFORCE);
1779 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All,
1780 forceOutNonbonded->forceWithShiftForces().force());
1781 wallcycle_start_nocount(wcycle, ewcFORCE);
1784 /* If there are multiple fshift output buffers we need to reduce them */
1785 if (stepWork.computeVirial)
1787 /* This is not in a subcounter because it takes a
1788 negligible and constant-sized amount of time */
1789 nbnxn_atomdata_add_nbat_fshift_to_fshift(
1790 *nbv->nbat, forceOutNonbonded->forceWithShiftForces().shiftForces());
1794 // TODO Force flags should include haveFreeEnergyWork for this domain
1795 if (stepWork.useGpuXHalo && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1797 wallcycle_stop(wcycle, ewcFORCE);
1798 /* Wait for non-local coordinate data to be copied from device */
1799 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1800 wallcycle_start_nocount(wcycle, ewcFORCE);
1803 // Compute wall interactions, when present.
1804 // Note: should be moved to special forces.
1805 if (inputrec->nwall && stepWork.computeNonbondedForces)
1807 /* foreign lambda component for walls */
1808 real dvdl_walls = do_walls(*inputrec,
1812 x.unpaddedConstArrayRef(),
1813 &forceOutMtsLevel0.forceWithVirial(),
1815 enerd->grpp.ener[egLJSR].data(),
1817 enerd->dvdl_lin[efptVDW] += dvdl_walls;
1820 if (stepWork.computeListedForces)
1822 /* Check whether we need to take into account PBC in listed interactions */
1823 bool needMolPbc = false;
1824 for (const auto& listedForces : fr->listedForces)
1826 if (listedForces.haveCpuListedForces(*fr->fcdata))
1828 needMolPbc = fr->bMolPBC;
1836 /* Since all atoms are in the rectangular or triclinic unit-cell,
1837 * only single box vector shifts (2 in x) are required.
1839 set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1842 for (int mtsIndex = 0; mtsIndex < (fr->useMts && stepWork.computeSlowForces ? 2 : 1); mtsIndex++)
1844 ListedForces& listedForces = fr->listedForces[mtsIndex];
1845 ForceOutputs& forceOut = (mtsIndex == 0 ? forceOutMtsLevel0 : *forceOutMtsLevel1);
1846 listedForces.calculate(wcycle,
1862 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr,
1867 if (stepWork.computeSlowForces)
1869 calculateLongRangeNonbondeds(fr,
1875 x.unpaddedConstArrayRef(),
1876 &forceOutMtsLevel1->forceWithVirial(),
1880 as_rvec_array(dipoleData.muStateAB),
1882 ddBalanceRegionHandler);
1885 wallcycle_stop(wcycle, ewcFORCE);
1887 // VdW dispersion correction, only computed on master rank to avoid double counting
1888 if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1890 // Calculate long range corrections to pressure and energy
1891 const DispersionCorrection::Correction correction =
1892 fr->dispersionCorrection->calculate(box, lambda[efptVDW]);
1894 if (stepWork.computeEnergy)
1896 enerd->term[F_DISPCORR] = correction.energy;
1897 enerd->term[F_DVDL_VDW] += correction.dvdl;
1898 enerd->dvdl_lin[efptVDW] += correction.dvdl;
1900 if (stepWork.computeVirial)
1902 correction.correctVirial(vir_force);
1903 enerd->term[F_PDISPCORR] = correction.pressure;
1907 computeSpecialForces(fplog,
1919 x.unpaddedArrayRef(),
1923 &forceOutMtsLevel0.forceWithVirial(),
1924 forceOutMtsLevel1 ? &forceOutMtsLevel1->forceWithVirial() : nullptr,
1927 stepWork.doNeighborSearch);
1929 if (havePPDomainDecomposition(cr) && stepWork.computeForces && stepWork.useGpuFHalo
1930 && domainWork.haveCpuLocalForceWork)
1932 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(), AtomLocality::Local);
1935 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
1936 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
1937 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFHalo),
1938 "The schedule below does not allow for nonbonded MTS with GPU halo exchange");
1939 // Will store the amount of cycles spent waiting for the GPU that
1940 // will be later used in the DLB accounting.
1941 float cycles_wait_gpu = 0;
1942 if (useOrEmulateGpuNb && stepWork.computeNonbondedForces)
1944 auto& forceWithShiftForces = forceOutNonbonded->forceWithShiftForces();
1946 /* wait for non-local forces (or calculate in emulation mode) */
1947 if (havePPDomainDecomposition(cr))
1949 if (simulationWork.useGpuNonbonded)
1951 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1953 AtomLocality::NonLocal,
1954 enerd->grpp.ener[egLJSR].data(),
1955 enerd->grpp.ener[egCOULSR].data(),
1956 forceWithShiftForces.shiftForces(),
1961 wallcycle_start_nocount(wcycle, ewcFORCE);
1963 fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes, step, nrnb, wcycle);
1964 wallcycle_stop(wcycle, ewcFORCE);
1967 if (stepWork.useGpuFBufferOps)
1969 // TODO: move this into DomainLifetimeWorkload, including the second part of the
1970 // condition The bonded and free energy CPU tasks can have non-local force
1971 // contributions which are a dependency for the GPU force reduction.
1972 bool haveNonLocalForceContribInCpuBuffer =
1973 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1975 if (haveNonLocalForceContribInCpuBuffer)
1977 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
1978 AtomLocality::NonLocal);
1982 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->execute();
1984 if (!stepWork.useGpuFHalo)
1986 // copy from GPU input for dd_move_f()
1987 stateGpu->copyForcesFromGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
1988 AtomLocality::NonLocal);
1993 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
1996 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
1998 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
2003 /* Combining the forces for multiple time stepping before the halo exchange, when possible,
2004 * avoids an extra halo exchange (when DD is used) and post-processing step.
2006 const bool combineMtsForcesBeforeHaloExchange =
2007 (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2008 && (legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0
2009 && !(stepWork.computeVirial || simulationWork.useGpuNonbonded || useGpuPmeOnThisRank));
2010 if (combineMtsForcesBeforeHaloExchange)
2012 const int numAtoms = havePPDomainDecomposition(cr) ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr;
2013 combineMtsForces(numAtoms,
2014 force.unpaddedArrayRef(),
2015 forceView->forceMtsCombined(),
2016 inputrec->mtsLevels[1].stepFactor);
2019 if (havePPDomainDecomposition(cr))
2021 /* We are done with the CPU compute.
2022 * We will now communicate the non-local forces.
2023 * If we use a GPU this will overlap with GPU work, so in that case
2024 * we do not close the DD force balancing region here.
2026 ddBalanceRegionHandler.closeAfterForceComputationCpu();
2028 if (stepWork.computeForces)
2031 if (stepWork.useGpuFHalo)
2033 communicateGpuHaloForces(*cr, domainWork.haveCpuLocalForceWork);
2037 if (stepWork.useGpuFBufferOps)
2039 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
2042 // Without MTS or with MTS at slow steps with uncombined forces we need to
2043 // communicate the fast forces
2044 if (!fr->useMts || !combineMtsForcesBeforeHaloExchange)
2046 dd_move_f(cr->dd, &forceOutMtsLevel0.forceWithShiftForces(), wcycle);
2048 // With MTS we need to communicate the slow or combined (in forceOutMtsLevel1) forces
2049 if (fr->useMts && stepWork.computeSlowForces)
2051 dd_move_f(cr->dd, &forceOutMtsLevel1->forceWithShiftForces(), wcycle);
2057 // With both nonbonded and PME offloaded a GPU on the same rank, we use
2058 // an alternating wait/reduction scheme.
2059 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
2060 && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
2061 if (alternateGpuWait)
2063 alternatePmeNbGpuWaitReduce(
2064 fr->nbv.get(), fr->pmedata, forceOutNonbonded, forceOutMtsLevel1, enerd, lambda[efptCOUL], stepWork, wcycle);
2067 if (!alternateGpuWait && useGpuPmeOnThisRank)
2069 pme_gpu_wait_and_reduce(
2070 fr->pmedata, stepWork, wcycle, &forceOutMtsLevel1->forceWithVirial(), enerd, lambda[efptCOUL]);
2073 /* Wait for local GPU NB outputs on the non-alternating wait path */
2074 if (!alternateGpuWait && stepWork.computeNonbondedForces && simulationWork.useGpuNonbonded)
2076 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
2077 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
2078 * but even with a step of 0.1 ms the difference is less than 1%
2081 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
2082 const float waitCycles =
2083 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
2085 AtomLocality::Local,
2086 enerd->grpp.ener[egLJSR].data(),
2087 enerd->grpp.ener[egCOULSR].data(),
2088 forceOutNonbonded->forceWithShiftForces().shiftForces(),
2091 if (ddBalanceRegionHandler.useBalancingRegion())
2093 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
2094 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
2096 /* We measured few cycles, it could be that the kernel
2097 * and transfer finished earlier and there was no actual
2098 * wait time, only API call overhead.
2099 * Then the actual time could be anywhere between 0 and
2100 * cycles_wait_est. We will use half of cycles_wait_est.
2102 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
2104 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
2108 if (fr->nbv->emulateGpu())
2110 // NOTE: emulation kernel is not included in the balancing region,
2111 // but emulation mode does not target performance anyway
2112 wallcycle_start_nocount(wcycle, ewcFORCE);
2117 InteractionLocality::Local,
2118 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
2122 wallcycle_stop(wcycle, ewcFORCE);
2125 // If on GPU PME-PP comms path, receive forces from PME before GPU buffer ops
2126 // TODO refactor this and unify with below default-path call to the same function
2127 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces
2128 && simulationWork.useGpuPmePpCommunication)
2130 /* In case of node-splitting, the PP nodes receive the long-range
2131 * forces, virial and energy from the PME nodes here.
2133 pme_receive_force_ener(fr,
2135 &forceOutMtsLevel1->forceWithVirial(),
2137 simulationWork.useGpuPmePpCommunication,
2138 stepWork.useGpuPmeFReduction,
2143 /* Do the nonbonded GPU (or emulation) force buffer reduction
2144 * on the non-alternating path. */
2145 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
2146 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
2147 if (useOrEmulateGpuNb && !alternateGpuWait)
2149 if (stepWork.useGpuFBufferOps)
2151 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2153 // Flag to specify whether the CPU force buffer has contributions to
2154 // local atoms. This depends on whether there are CPU-based force tasks
2155 // or when DD is active the halo exchange has resulted in contributions
2156 // from the non-local part.
2157 const bool haveLocalForceContribInCpuBuffer =
2158 (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
2160 // TODO: move these steps as early as possible:
2161 // - CPU f H2D should be as soon as all CPU-side forces are done
2162 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
2163 // before the next CPU task that consumes the forces: vsite spread or update)
2164 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
2165 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
2166 // These should be unified.
2167 if (haveLocalForceContribInCpuBuffer && !stepWork.useGpuFHalo)
2169 // Note: AtomLocality::All is used for the non-DD case because, as in this
2170 // case copyForcesToGpu() uses a separate stream, it allows overlap of
2171 // CPU force H2D with GPU force tasks on all streams including those in the
2172 // local stream which would otherwise be implicit dependencies for the
2173 // transfer and would not overlap.
2174 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
2176 stateGpu->copyForcesToGpu(forceWithShift, locality);
2179 if (stepWork.computeNonbondedForces)
2181 fr->gpuForceReduction[gmx::AtomLocality::Local]->execute();
2184 // Copy forces to host if they are needed for update or if virtual sites are enabled.
2185 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
2186 // TODO: When the output flags will be included in step workload, this copy can be combined with the
2187 // copy call done in sim_utils(...) for the output.
2188 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
2189 // they should not be copied in do_md(...) for the output.
2190 if (!simulationWork.useGpuUpdate
2191 || (simulationWork.useGpuUpdate && DOMAINDECOMP(cr) && haveHostPmePpComms) || vsite)
2193 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
2194 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
2197 else if (stepWork.computeNonbondedForces)
2199 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2200 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
2204 launchGpuEndOfStepTasks(
2205 nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork, useGpuPmeOnThisRank, step, wcycle);
2207 if (DOMAINDECOMP(cr))
2209 dd_force_flop_stop(cr->dd, nrnb);
2212 const bool haveCombinedMtsForces = (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2213 && combineMtsForcesBeforeHaloExchange);
2214 if (stepWork.computeForces)
2216 postProcessForceWithShiftForces(
2217 nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutMtsLevel0, vir_force, *mdatoms, *fr, vsite, stepWork);
2219 if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2221 postProcessForceWithShiftForces(
2222 nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, *mdatoms, *fr, vsite, stepWork);
2226 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
2227 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
2228 && stepWork.computeSlowForces)
2230 /* In case of node-splitting, the PP nodes receive the long-range
2231 * forces, virial and energy from the PME nodes here.
2233 pme_receive_force_ener(fr,
2235 &forceOutMtsLevel1->forceWithVirial(),
2237 simulationWork.useGpuPmePpCommunication,
2242 if (stepWork.computeForces)
2244 /* If we don't use MTS or if we already combined the MTS forces before, we only
2245 * need to post-process one ForceOutputs object here, called forceOutCombined,
2246 * otherwise we have to post-process two outputs and then combine them.
2248 ForceOutputs& forceOutCombined = (haveCombinedMtsForces ? forceOutMts.value() : forceOutMtsLevel0);
2250 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutCombined, vir_force, mdatoms, fr, vsite, stepWork);
2252 if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2255 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, mdatoms, fr, vsite, stepWork);
2257 combineMtsForces(mdatoms->homenr,
2258 force.unpaddedArrayRef(),
2259 forceView->forceMtsCombined(),
2260 inputrec->mtsLevels[1].stepFactor);
2264 if (stepWork.computeEnergy)
2266 /* Compute the final potential energy terms */
2267 accumulatePotentialEnergies(enerd, lambda, inputrec->fepvals);
2269 if (!EI_TPI(inputrec->eI))
2271 checkPotentialEnergyValidity(step, *enerd, *inputrec);
2275 /* In case we don't have constraints and are using GPUs, the next balancing
2276 * region starts here.
2277 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
2278 * virial calculation and COM pulling, is not thus not included in
2279 * the balance timing, which is ok as most tasks do communication.
2281 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);