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/nonbonded.h"
62 #include "gromacs/gmxlib/nrnb.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,
208 gmx::arrayRefFromArray(mdatoms->massT, mdatoms->nr),
212 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Restraint)],
216 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Restraint] += dvdl;
217 wallcycle_stop(wcycle, ewcPULLPOT);
220 static void pme_receive_force_ener(t_forcerec* fr,
222 gmx::ForceWithVirial* forceWithVirial,
223 gmx_enerdata_t* enerd,
224 bool useGpuPmePpComms,
225 bool receivePmeForceToGpu,
226 gmx_wallcycle_t wcycle)
228 real e_q, e_lj, dvdl_q, dvdl_lj;
229 float cycles_ppdpme, cycles_seppme;
231 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
232 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
234 /* In case of node-splitting, the PP nodes receive the long-range
235 * forces, virial and energy from the PME nodes here.
237 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
240 gmx_pme_receive_f(fr->pmePpCommGpu.get(),
248 receivePmeForceToGpu,
250 enerd->term[F_COUL_RECIP] += e_q;
251 enerd->term[F_LJ_RECIP] += e_lj;
252 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Coul] += dvdl_q;
253 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += dvdl_lj;
257 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
259 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
262 static void print_large_forces(FILE* fp,
267 ArrayRef<const RVec> x,
268 ArrayRef<const RVec> f)
270 real force2Tolerance = gmx::square(forceTolerance);
271 gmx::index numNonFinite = 0;
272 for (int i = 0; i < md->homenr; i++)
274 real force2 = norm2(f[i]);
275 bool nonFinite = !std::isfinite(force2);
276 if (force2 >= force2Tolerance || nonFinite)
279 "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
292 if (numNonFinite > 0)
294 /* Note that with MPI this fatal call on one rank might interrupt
295 * the printing on other ranks. But we can only avoid that with
296 * an expensive MPI barrier that we would need at each step.
298 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
302 //! When necessary, spreads forces on vsites and computes the virial for \p forceOutputs->forceWithShiftForces()
303 static void postProcessForceWithShiftForces(t_nrnb* nrnb,
304 gmx_wallcycle_t wcycle,
306 ArrayRef<const RVec> x,
307 ForceOutputs* forceOutputs,
309 const t_mdatoms& mdatoms,
310 const t_forcerec& fr,
311 gmx::VirtualSitesHandler* vsite,
312 const StepWorkload& stepWork)
314 ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
316 /* If we have NoVirSum forces, but we do not calculate the virial,
317 * we later sum the forceWithShiftForces buffer together with
318 * the noVirSum buffer and spread the combined vsite forces at once.
320 if (vsite && (!forceOutputs->haveForceWithVirial() || stepWork.computeVirial))
322 using VirialHandling = gmx::VirtualSitesHandler::VirialHandling;
324 auto f = forceWithShiftForces.force();
325 auto fshift = forceWithShiftForces.shiftForces();
326 const VirialHandling virialHandling =
327 (stepWork.computeVirial ? VirialHandling::Pbc : VirialHandling::None);
328 vsite->spreadForces(x, f, virialHandling, fshift, nullptr, nrnb, box, wcycle);
329 forceWithShiftForces.haveSpreadVsiteForces() = true;
332 if (stepWork.computeVirial)
334 /* Calculation of the virial must be done after vsites! */
336 0, mdatoms.homenr, as_rvec_array(x.data()), forceWithShiftForces, vir_force, box, nrnb, &fr, fr.pbcType);
340 //! Spread, compute virial for and sum forces, when necessary
341 static void postProcessForces(const t_commrec* cr,
344 gmx_wallcycle_t wcycle,
346 ArrayRef<const RVec> x,
347 ForceOutputs* forceOutputs,
349 const t_mdatoms* mdatoms,
350 const t_forcerec* fr,
351 gmx::VirtualSitesHandler* vsite,
352 const StepWorkload& stepWork)
354 // Extract the final output force buffer, which is also the buffer for forces with shift forces
355 ArrayRef<RVec> f = forceOutputs->forceWithShiftForces().force();
357 if (forceOutputs->haveForceWithVirial())
359 auto& forceWithVirial = forceOutputs->forceWithVirial();
363 /* Spread the mesh force on virtual sites to the other particles...
364 * This is parallellized. MPI communication is performed
365 * if the constructing atoms aren't local.
367 GMX_ASSERT(!stepWork.computeVirial || f.data() != forceWithVirial.force_.data(),
368 "We need separate force buffers for shift and virial forces when "
369 "computing the virial");
370 GMX_ASSERT(!stepWork.computeVirial
371 || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
372 "We should spread the force with shift forces separately when computing "
374 const gmx::VirtualSitesHandler::VirialHandling virialHandling =
375 (stepWork.computeVirial ? gmx::VirtualSitesHandler::VirialHandling::NonLinear
376 : gmx::VirtualSitesHandler::VirialHandling::None);
377 matrix virial = { { 0 } };
378 vsite->spreadForces(x, forceWithVirial.force_, virialHandling, {}, virial, nrnb, box, wcycle);
379 forceWithVirial.addVirialContribution(virial);
382 if (stepWork.computeVirial)
384 /* Now add the forces, this is local */
385 sum_forces(f, forceWithVirial.force_);
387 /* Add the direct virial contributions */
389 forceWithVirial.computeVirial_,
390 "forceWithVirial should request virial computation when we request the virial");
391 m_add(vir_force, forceWithVirial.getVirial(), vir_force);
395 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
401 GMX_ASSERT(vsite == nullptr || forceOutputs->forceWithShiftForces().haveSpreadVsiteForces(),
402 "We should have spread the vsite forces (earlier)");
405 if (fr->print_force >= 0)
407 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
411 static void do_nb_verlet(t_forcerec* fr,
412 const interaction_const_t* ic,
413 gmx_enerdata_t* enerd,
414 const StepWorkload& stepWork,
415 const InteractionLocality ilocality,
419 gmx_wallcycle_t wcycle)
421 if (!stepWork.computeNonbondedForces)
423 /* skip non-bonded calculation */
427 nonbonded_verlet_t* nbv = fr->nbv.get();
429 /* GPU kernel launch overhead is already timed separately */
432 /* When dynamic pair-list pruning is requested, we need to prune
433 * at nstlistPrune steps.
435 if (nbv->isDynamicPruningStepCpu(step))
437 /* Prune the pair-list beyond fr->ic->rlistPrune using
438 * the current coordinates of the atoms.
440 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
441 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
442 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
446 nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
449 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
451 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, v.ssize());
453 /* Note that we would like to avoid this conditional by putting it
454 * into the omp pragma instead, but then we still take the full
455 * omp parallel for overhead (at least with gcc5).
457 if (!useOpenmpThreading || nth == 1)
466 #pragma omp parallel for num_threads(nth) schedule(static)
467 for (gmx::index i = 0; i < v.ssize(); i++)
474 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
476 * \param groupOptions Group options, containing T-coupling options
478 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
480 real nrdfCoupled = 0;
481 real nrdfUncoupled = 0;
482 real kineticEnergy = 0;
483 for (int g = 0; g < groupOptions.ngtc; g++)
485 if (groupOptions.tau_t[g] >= 0)
487 nrdfCoupled += groupOptions.nrdf[g];
488 kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * gmx::c_boltz;
492 nrdfUncoupled += groupOptions.nrdf[g];
496 /* This conditional with > also catches nrdf=0 */
497 if (nrdfCoupled > nrdfUncoupled)
499 return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
507 /*! \brief This routine checks that the potential energy is finite.
509 * Always checks that the potential energy is finite. If step equals
510 * inputrec.init_step also checks that the magnitude of the potential energy
511 * is reasonable. Terminates with a fatal error when a check fails.
512 * Note that passing this check does not guarantee finite forces,
513 * since those use slightly different arithmetics. But in most cases
514 * there is just a narrow coordinate range where forces are not finite
515 * and energies are finite.
517 * \param[in] step The step number, used for checking and printing
518 * \param[in] enerd The energy data; the non-bonded group energies need to be added to
519 * enerd.term[F_EPOT] before calling this routine \param[in] inputrec The input record
521 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
523 /* Threshold valid for comparing absolute potential energy against
524 * the kinetic energy. Normally one should not consider absolute
525 * potential energy values, but with a factor of one million
526 * we should never get false positives.
528 constexpr real c_thresholdFactor = 1e6;
530 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
531 real averageKineticEnergy = 0;
532 /* We only check for large potential energy at the initial step,
533 * because that is by far the most likely step for this too occur
534 * and because computing the average kinetic energy is not free.
535 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
536 * before they become NaN.
538 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
540 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
543 if (energyIsNotFinite
544 || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
549 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
550 "contributions to the energy are %g and %g, respectively. A %s potential energy "
551 "can be caused by overlapping interactions in bonded interactions or very large%s "
552 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
553 "configuration, incorrect interactions or parameters in the topology.",
556 energyIsNotFinite ? "not finite" : "extremely high",
558 enerd.term[F_COUL_SR],
559 energyIsNotFinite ? "non-finite" : "very high",
560 energyIsNotFinite ? " or Nan" : "");
564 /*! \brief Return true if there are special forces computed this step.
566 * The conditionals exactly correspond to those in computeSpecialForces().
568 static bool haveSpecialForces(const t_inputrec& inputrec,
569 const gmx::ForceProviders& forceProviders,
570 const pull_t* pull_work,
571 const bool computeForces,
575 return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
576 (inputrec.bPull && pull_have_potential(*pull_work)) || // pull
577 inputrec.bRot || // enforced rotation
578 (ed != nullptr) || // flooding
579 (inputrec.bIMD && computeForces)); // IMD
582 /*! \brief Compute forces and/or energies for special algorithms
584 * The intention is to collect all calls to algorithms that compute
585 * forces on local atoms only and that do not contribute to the local
586 * virial sum (but add their virial contribution separately).
587 * Eventually these should likely all become ForceProviders.
588 * Within this function the intention is to have algorithms that do
589 * global communication at the end, so global barriers within the MD loop
590 * are as close together as possible.
592 * \param[in] fplog The log file
593 * \param[in] cr The communication record
594 * \param[in] inputrec The input record
595 * \param[in] awh The Awh module (nullptr if none in use).
596 * \param[in] enforcedRotation Enforced rotation module.
597 * \param[in] imdSession The IMD session
598 * \param[in] pull_work The pull work structure.
599 * \param[in] step The current MD step
600 * \param[in] t The current time
601 * \param[in,out] wcycle Wallcycle accounting struct
602 * \param[in,out] forceProviders Pointer to a list of force providers
603 * \param[in] box The unit cell
604 * \param[in] x The coordinates
605 * \param[in] mdatoms Per atom properties
606 * \param[in] lambda Array of free-energy lambda values
607 * \param[in] stepWork Step schedule flags
608 * \param[in,out] forceWithVirialMtsLevel0 Force and virial for MTS level0 forces
609 * \param[in,out] forceWithVirialMtsLevel1 Force and virial for MTS level1 forces, can be nullptr
610 * \param[in,out] enerd Energy buffer
611 * \param[in,out] ed Essential dynamics pointer
612 * \param[in] didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
614 * \todo Remove didNeighborSearch, which is used incorrectly.
615 * \todo Convert all other algorithms called here to ForceProviders.
617 static void computeSpecialForces(FILE* fplog,
619 const t_inputrec& inputrec,
621 gmx_enfrot* enforcedRotation,
622 gmx::ImdSession* imdSession,
626 gmx_wallcycle_t wcycle,
627 gmx::ForceProviders* forceProviders,
629 gmx::ArrayRef<const gmx::RVec> x,
630 const t_mdatoms* mdatoms,
631 gmx::ArrayRef<const real> lambda,
632 const StepWorkload& stepWork,
633 gmx::ForceWithVirial* forceWithVirialMtsLevel0,
634 gmx::ForceWithVirial* forceWithVirialMtsLevel1,
635 gmx_enerdata_t* enerd,
637 bool didNeighborSearch)
639 /* NOTE: Currently all ForceProviders only provide forces.
640 * When they also provide energies, remove this conditional.
642 if (stepWork.computeForces)
644 gmx::ForceProviderInput forceProviderInput(
647 gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->homenr),
648 gmx::arrayRefFromArray(mdatoms->massT, mdatoms->homenr),
652 gmx::ForceProviderOutput forceProviderOutput(forceWithVirialMtsLevel0, enerd);
654 /* Collect forces from modules */
655 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
658 if (inputrec.bPull && pull_have_potential(*pull_work))
660 const int mtsLevel = forceGroupMtsLevel(inputrec.mtsLevels, gmx::MtsForceGroups::Pull);
661 if (mtsLevel == 0 || stepWork.computeSlowForces)
663 auto& forceWithVirial = (mtsLevel == 0) ? forceWithVirialMtsLevel0 : forceWithVirialMtsLevel1;
664 pull_potential_wrapper(
665 cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work, lambda.data(), t, wcycle);
670 const int mtsLevel = forceGroupMtsLevel(inputrec.mtsLevels, gmx::MtsForceGroups::Pull);
671 if (mtsLevel == 0 || stepWork.computeSlowForces)
673 const bool needForeignEnergyDifferences = awh->needForeignEnergyDifferences(step);
674 std::vector<double> foreignLambdaDeltaH, foreignLambdaDhDl;
675 if (needForeignEnergyDifferences)
677 enerd->foreignLambdaTerms.finalizePotentialContributions(
678 enerd->dvdl_lin, lambda, *inputrec.fepvals);
679 std::tie(foreignLambdaDeltaH, foreignLambdaDhDl) = enerd->foreignLambdaTerms.getTerms(cr);
682 auto& forceWithVirial = (mtsLevel == 0) ? forceWithVirialMtsLevel0 : forceWithVirialMtsLevel1;
683 enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(
685 gmx::arrayRefFromArray(mdatoms->massT, mdatoms->nr),
696 /* Add the forces from enforced rotation potentials (if any) */
699 wallcycle_start(wcycle, ewcROTadd);
700 enerd->term[F_COM_PULL] +=
701 add_rot_forces(enforcedRotation, forceWithVirialMtsLevel0->force_, cr, step, t);
702 wallcycle_stop(wcycle, ewcROTadd);
707 /* Note that since init_edsam() is called after the initialization
708 * of forcerec, edsam doesn't request the noVirSum force buffer.
709 * Thus if no other algorithm (e.g. PME) requires it, the forces
710 * here will contribute to the virial.
712 do_flood(cr, inputrec, x, forceWithVirialMtsLevel0->force_, ed, box, step, didNeighborSearch);
715 /* Add forces from interactive molecular dynamics (IMD), if any */
716 if (inputrec.bIMD && stepWork.computeForces)
718 imdSession->applyForces(forceWithVirialMtsLevel0->force_);
722 /*! \brief Launch the prepare_step and spread stages of PME GPU.
724 * \param[in] pmedata The PME structure
725 * \param[in] box The box matrix
726 * \param[in] stepWork Step schedule flags
727 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in the device memory.
728 * \param[in] lambdaQ The Coulomb lambda of the current state.
729 * \param[in] wcycle The wallcycle structure
731 static inline void launchPmeGpuSpread(gmx_pme_t* pmedata,
733 const StepWorkload& stepWork,
734 GpuEventSynchronizer* xReadyOnDevice,
736 gmx_wallcycle_t wcycle)
738 pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
739 pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle, lambdaQ);
742 /*! \brief Launch the FFT and gather stages of PME GPU
744 * This function only implements setting the output forces (no accumulation).
746 * \param[in] pmedata The PME structure
747 * \param[in] lambdaQ The Coulomb lambda of the current system state.
748 * \param[in] wcycle The wallcycle structure
749 * \param[in] stepWork Step schedule flags
751 static void launchPmeGpuFftAndGather(gmx_pme_t* pmedata,
753 gmx_wallcycle_t wcycle,
754 const gmx::StepWorkload& stepWork)
756 pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
757 pme_gpu_launch_gather(pmedata, wcycle, lambdaQ);
761 * Polling wait for either of the PME or nonbonded GPU tasks.
763 * Instead of a static order in waiting for GPU tasks, this function
764 * polls checking which of the two tasks completes first, and does the
765 * associated force buffer reduction overlapped with the other task.
766 * By doing that, unlike static scheduling order, it can always overlap
767 * one of the reductions, regardless of the GPU task completion order.
769 * \param[in] nbv Nonbonded verlet structure
770 * \param[in,out] pmedata PME module data
771 * \param[in,out] forceOutputsNonbonded Force outputs for the non-bonded forces and shift forces
772 * \param[in,out] forceOutputsPme Force outputs for the PME forces and virial
773 * \param[in,out] enerd Energy data structure results are reduced into
774 * \param[in] lambdaQ The Coulomb lambda of the current system state.
775 * \param[in] stepWork Step schedule flags
776 * \param[in] wcycle The wallcycle structure
778 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
780 gmx::ForceOutputs* forceOutputsNonbonded,
781 gmx::ForceOutputs* forceOutputsPme,
782 gmx_enerdata_t* enerd,
784 const StepWorkload& stepWork,
785 gmx_wallcycle_t wcycle)
787 bool isPmeGpuDone = false;
788 bool isNbGpuDone = false;
790 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
792 while (!isPmeGpuDone || !isNbGpuDone)
796 GpuTaskCompletion completionType =
797 (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
798 isPmeGpuDone = pme_gpu_try_finish_task(
799 pmedata, stepWork, wcycle, &forceOutputsPme->forceWithVirial(), enerd, lambdaQ, completionType);
804 auto& forceBuffersNonbonded = forceOutputsNonbonded->forceWithShiftForces();
805 GpuTaskCompletion completionType =
806 (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
807 isNbGpuDone = Nbnxm::gpu_try_finish_task(
811 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
812 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
813 forceBuffersNonbonded.shiftForces(),
819 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceBuffersNonbonded.force());
825 /*! \brief Set up the different force buffers; also does clearing.
827 * \param[in] forceHelperBuffers Helper force buffers
828 * \param[in] force force array
829 * \param[in] domainWork Domain lifetime workload flags
830 * \param[in] stepWork Step schedule flags
831 * \param[in] havePpDomainDecomposition Whether we have a PP domain decomposition
832 * \param[out] wcycle wallcycle recording structure
834 * \returns Cleared force output structure
836 static ForceOutputs setupForceOutputs(ForceHelperBuffers* forceHelperBuffers,
837 gmx::ArrayRefWithPadding<gmx::RVec> force,
838 const DomainLifetimeWorkload& domainWork,
839 const StepWorkload& stepWork,
840 const bool havePpDomainDecomposition,
841 gmx_wallcycle_t wcycle)
843 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
845 /* NOTE: We assume fr->shiftForces is all zeros here */
846 gmx::ForceWithShiftForces forceWithShiftForces(
847 force, stepWork.computeVirial, forceHelperBuffers->shiftForces());
849 if (stepWork.computeForces
850 && (domainWork.haveCpuLocalForceWork || !stepWork.useGpuFBufferOps
851 || (havePpDomainDecomposition && !stepWork.useGpuFHalo)))
853 /* Clear the short- and long-range forces */
854 clearRVecs(forceWithShiftForces.force(), true);
856 /* Clear the shift forces */
857 clearRVecs(forceWithShiftForces.shiftForces(), false);
860 /* If we need to compute the virial, we might need a separate
861 * force buffer for algorithms for which the virial is calculated
862 * directly, such as PME. Otherwise, forceWithVirial uses the
863 * the same force (f in legacy calls) buffer as other algorithms.
865 const bool useSeparateForceWithVirialBuffer =
866 (stepWork.computeForces
867 && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
868 /* forceWithVirial uses the local atom range only */
869 gmx::ForceWithVirial forceWithVirial(
870 useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
871 : force.unpaddedArrayRef(),
872 stepWork.computeVirial);
874 if (useSeparateForceWithVirialBuffer)
876 /* TODO: update comment
877 * We only compute forces on local atoms. Note that vsites can
878 * spread to non-local atoms, but that part of the buffer is
879 * cleared separately in the vsite spreading code.
881 clearRVecs(forceWithVirial.force_, true);
884 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
887 forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(), forceWithVirial);
891 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
893 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec& inputrec,
894 const t_forcerec& fr,
895 const pull_t* pull_work,
897 const t_mdatoms& mdatoms,
898 const SimulationWorkload& simulationWork,
899 const StepWorkload& stepWork)
901 DomainLifetimeWorkload domainWork;
902 // Note that haveSpecialForces is constant over the whole run
903 domainWork.haveSpecialForces =
904 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
905 domainWork.haveCpuListedForceWork = false;
906 domainWork.haveCpuBondedWork = false;
907 for (const auto& listedForces : fr.listedForces)
909 if (listedForces.haveCpuListedForces(*fr.fcdata))
911 domainWork.haveCpuListedForceWork = true;
913 if (listedForces.haveCpuBondeds())
915 domainWork.haveCpuBondedWork = true;
918 domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
919 // Note that haveFreeEnergyWork is constant over the whole run
920 domainWork.haveFreeEnergyWork =
921 (fr.efep != FreeEnergyPerturbationType::No && mdatoms.nPerturbed != 0);
922 // We assume we have local force work if there are CPU
923 // force tasks including PME or nonbondeds.
924 domainWork.haveCpuLocalForceWork =
925 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
926 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
927 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
932 /*! \brief Set up force flag stuct from the force bitmask.
934 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
935 * \param[in] mtsLevels The multiple time-stepping levels, either empty or 2 levels
936 * \param[in] step The current MD step
937 * \param[in] simulationWork Simulation workload description.
938 * \param[in] rankHasPmeDuty If this rank computes PME.
940 * \returns New Stepworkload description.
942 static StepWorkload setupStepWorkload(const int legacyFlags,
943 ArrayRef<const gmx::MtsLevel> mtsLevels,
945 const SimulationWorkload& simulationWork,
946 const bool rankHasPmeDuty)
948 GMX_ASSERT(mtsLevels.empty() || mtsLevels.size() == 2, "Expect 0 or 2 MTS levels");
949 const bool computeSlowForces = (mtsLevels.empty() || step % mtsLevels[1].stepFactor == 0);
952 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
953 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
954 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
955 flags.computeSlowForces = computeSlowForces;
956 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
957 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
958 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
959 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
960 flags.computeNonbondedForces =
961 ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded
962 && !(simulationWork.computeNonbondedAtMtsLevel1 && !computeSlowForces);
963 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
965 if (simulationWork.useGpuBufferOps)
967 GMX_ASSERT(simulationWork.useGpuNonbonded,
968 "Can only offload buffer ops if nonbonded computation is also offloaded");
970 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
971 // on virial steps the CPU reduction path is taken
972 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
973 flags.useGpuPmeFReduction = flags.computeSlowForces && flags.useGpuFBufferOps && simulationWork.useGpuPme
974 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication);
975 flags.useGpuXHalo = simulationWork.useGpuHaloExchange;
976 flags.useGpuFHalo = simulationWork.useGpuHaloExchange && flags.useGpuFBufferOps;
982 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
984 * TODO: eliminate \p useGpuPmeOnThisRank when this is
985 * incorporated in DomainLifetimeWorkload.
987 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
988 gmx::GpuBonded* gpuBonded,
990 gmx_enerdata_t* enerd,
991 const gmx::MdrunScheduleWorkload& runScheduleWork,
992 bool useGpuPmeOnThisRank,
994 gmx_wallcycle_t wcycle)
996 if (runScheduleWork.simulationWork.useGpuNonbonded && runScheduleWork.stepWork.computeNonbondedForces)
998 /* Launch pruning before buffer clearing because the API overhead of the
999 * clear kernel launches can leave the GPU idle while it could be running
1002 if (nbv->isDynamicPruningStepGpu(step))
1004 nbv->dispatchPruneKernelGpu(step);
1007 /* now clear the GPU outputs while we finish the step on the CPU */
1008 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1009 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1010 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
1011 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1012 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1015 if (useGpuPmeOnThisRank)
1017 pme_gpu_reinit_computation(pmedata, wcycle);
1020 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
1022 // in principle this should be included in the DD balancing region,
1023 // but generally it is infrequent so we'll omit it for the sake of
1025 gpuBonded->waitAccumulateEnergyTerms(enerd);
1027 gpuBonded->clearEnergies();
1031 //! \brief Data structure to hold dipole-related data and staging arrays
1034 //! Dipole staging for fast summing over MPI
1035 gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
1036 //! Dipole staging for states A and B (index 0 and 1 resp.)
1037 gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
1041 static void reduceAndUpdateMuTot(DipoleData* dipoleData,
1042 const t_commrec* cr,
1043 const bool haveFreeEnergy,
1044 gmx::ArrayRef<const real> lambda,
1046 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1050 gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
1051 ddBalanceRegionHandler.reopenRegionCpu();
1053 for (int i = 0; i < 2; i++)
1055 for (int j = 0; j < DIM; j++)
1057 dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
1061 if (!haveFreeEnergy)
1063 copy_rvec(dipoleData->muStateAB[0], muTotal);
1067 for (int j = 0; j < DIM; j++)
1069 muTotal[j] = (1.0 - lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)])
1070 * dipoleData->muStateAB[0][j]
1071 + lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]
1072 * dipoleData->muStateAB[1][j];
1077 /*! \brief Combines MTS level0 and level1 force buffes into a full and MTS-combined force buffer.
1079 * \param[in] numAtoms The number of atoms to combine forces for
1080 * \param[in,out] forceMtsLevel0 Input: F_level0, output: F_level0 + F_level1
1081 * \param[in,out] forceMts Input: F_level1, output: F_level0 + mtsFactor * F_level1
1082 * \param[in] mtsFactor The factor between the level0 and level1 time step
1084 static void combineMtsForces(const int numAtoms,
1085 ArrayRef<RVec> forceMtsLevel0,
1086 ArrayRef<RVec> forceMts,
1087 const real mtsFactor)
1089 const int gmx_unused numThreads = gmx_omp_nthreads_get(emntDefault);
1090 #pragma omp parallel for num_threads(numThreads) schedule(static)
1091 for (int i = 0; i < numAtoms; i++)
1093 const RVec forceMtsLevel0Tmp = forceMtsLevel0[i];
1094 forceMtsLevel0[i] += forceMts[i];
1095 forceMts[i] = forceMtsLevel0Tmp + mtsFactor * forceMts[i];
1099 /*! \brief Setup for the local and non-local GPU force reductions:
1100 * reinitialization plus the registration of forces and dependencies.
1102 * \param [in] runScheduleWork Schedule workload flag structure
1103 * \param [in] cr Communication record object
1104 * \param [in] fr Force record object
1106 static void setupGpuForceReductions(gmx::MdrunScheduleWorkload* runScheduleWork,
1107 const t_commrec* cr,
1111 nonbonded_verlet_t* nbv = fr->nbv.get();
1112 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1114 // (re-)initialize local GPU force reduction
1115 const bool accumulate =
1116 runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr);
1117 const int atomStart = 0;
1118 fr->gpuForceReduction[gmx::AtomLocality::Local]->reinit(stateGpu->getForces(),
1119 nbv->getNumAtoms(AtomLocality::Local),
1120 nbv->getGridIndices(),
1123 stateGpu->fReducedOnDevice());
1125 // register forces and add dependencies
1126 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerNbnxmForce(nbv->getGpuForces());
1128 if (runScheduleWork->simulationWork.useGpuPme
1129 && (thisRankHasDuty(cr, DUTY_PME) || runScheduleWork->simulationWork.useGpuPmePpCommunication))
1131 DeviceBuffer<gmx::RVec> forcePtr =
1132 thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1133 : // PME force buffer on same GPU
1134 fr->pmePpCommGpu->getGpuForceStagingPtr(); // buffer received from other GPU
1135 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerRvecForce(forcePtr);
1137 GpuEventSynchronizer* const pmeSynchronizer =
1138 (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1139 : // PME force buffer on same GPU
1140 fr->pmePpCommGpu->getForcesReadySynchronizer()); // buffer received from other GPU
1141 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(pmeSynchronizer);
1144 if ((runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr))
1145 && !runScheduleWork->simulationWork.useGpuHaloExchange)
1147 auto forcesReadyLocality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
1148 const bool useGpuForceBufferOps = true;
1149 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1150 stateGpu->getForcesReadyOnDeviceEvent(forcesReadyLocality, useGpuForceBufferOps));
1153 if (runScheduleWork->simulationWork.useGpuHaloExchange)
1155 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1156 cr->dd->gpuHaloExchange[0][0]->getForcesReadyOnDeviceEvent());
1159 if (havePPDomainDecomposition(cr))
1161 // (re-)initialize non-local GPU force reduction
1162 const bool accumulate = runScheduleWork->domainWork.haveCpuBondedWork
1163 || runScheduleWork->domainWork.haveFreeEnergyWork;
1164 const int atomStart = dd_numHomeAtoms(*cr->dd);
1165 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->reinit(stateGpu->getForces(),
1166 nbv->getNumAtoms(AtomLocality::NonLocal),
1167 nbv->getGridIndices(),
1171 // register forces and add dependencies
1172 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->registerNbnxmForce(nbv->getGpuForces());
1173 if (runScheduleWork->domainWork.haveCpuBondedWork || runScheduleWork->domainWork.haveFreeEnergyWork)
1175 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->addDependency(
1176 stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::NonLocal, true));
1182 void do_force(FILE* fplog,
1183 const t_commrec* cr,
1184 const gmx_multisim_t* ms,
1185 const t_inputrec& inputrec,
1187 gmx_enfrot* enforcedRotation,
1188 gmx::ImdSession* imdSession,
1192 gmx_wallcycle_t wcycle,
1193 const gmx_localtop_t* top,
1195 gmx::ArrayRefWithPadding<gmx::RVec> x,
1196 const history_t* hist,
1197 gmx::ForceBuffersView* forceView,
1199 const t_mdatoms* mdatoms,
1200 gmx_enerdata_t* enerd,
1201 gmx::ArrayRef<const real> lambda,
1203 gmx::MdrunScheduleWorkload* runScheduleWork,
1204 gmx::VirtualSitesHandler* vsite,
1209 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1211 auto force = forceView->forceWithPadding();
1212 GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1213 "The size of the force buffer should be at least the number of atoms to compute "
1216 nonbonded_verlet_t* nbv = fr->nbv.get();
1217 interaction_const_t* ic = fr->ic.get();
1219 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1221 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1223 runScheduleWork->stepWork = setupStepWorkload(
1224 legacyFlags, inputrec.mtsLevels, step, simulationWork, thisRankHasDuty(cr, DUTY_PME));
1225 const StepWorkload& stepWork = runScheduleWork->stepWork;
1227 const bool useGpuPmeOnThisRank =
1228 simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces;
1230 /* At a search step we need to start the first balancing region
1231 * somewhere early inside the step after communication during domain
1232 * decomposition (and not during the previous step as usual).
1234 if (stepWork.doNeighborSearch)
1236 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1239 clear_mat(vir_force);
1241 if (fr->pbcType != PbcType::No)
1243 /* Compute shift vectors every step,
1244 * because of pressure coupling or box deformation!
1246 if (stepWork.haveDynamicBox && stepWork.stateChanged)
1248 calc_shifts(box, fr->shift_vec);
1251 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1252 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1255 put_atoms_in_box_omp(fr->pbcType,
1257 x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1258 gmx_omp_nthreads_get(emntDefault));
1259 inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1263 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1265 const bool pmeSendCoordinatesFromGpu =
1266 GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1267 const bool reinitGpuPmePpComms =
1268 GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1270 const auto localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1271 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1272 AtomLocality::Local, simulationWork, stepWork)
1275 // Copy coordinate from the GPU if update is on the GPU and there
1276 // are forces to be computed on the CPU, or for the computation of
1277 // virial, or if host-side data will be transferred from this task
1278 // to a remote task for halo exchange or PME-PP communication. At
1279 // search steps the current coordinates are already on the host,
1280 // hence copy is not needed.
1281 const bool haveHostPmePpComms =
1282 !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1284 GMX_ASSERT(simulationWork.useGpuHaloExchange
1285 == ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange[0].empty())),
1286 "The GPU halo exchange is active, but it has not been constructed.");
1287 const bool haveHostHaloExchangeComms =
1288 havePPDomainDecomposition(cr) && !simulationWork.useGpuHaloExchange;
1290 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1291 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1292 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1293 || haveHostPmePpComms || haveHostHaloExchangeComms))
1295 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1296 haveCopiedXFromGpu = true;
1299 // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1300 // Otherwise the send will occur after H2D coordinate transfer.
1301 if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu && stepWork.computeSlowForces)
1303 /* Send particle coordinates to the pme nodes */
1304 if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1306 GMX_ASSERT(haveCopiedXFromGpu,
1307 "a wait should only be triggered if copy has been scheduled");
1308 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1311 gmx_pme_send_coordinates(fr,
1314 as_rvec_array(x.unpaddedArrayRef().data()),
1315 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1316 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1317 (stepWork.computeVirial || stepWork.computeEnergy),
1319 simulationWork.useGpuPmePpCommunication,
1320 reinitGpuPmePpComms,
1321 pmeSendCoordinatesFromGpu,
1322 localXReadyOnDevice,
1326 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1327 // The local coordinates can be copied right away.
1328 // NOTE: Consider moving this copy to right after they are updated and constrained,
1329 // if the later is not offloaded.
1330 if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1332 if (stepWork.doNeighborSearch)
1334 // TODO refactor this to do_md, after partitioning.
1335 stateGpu->reinit(mdatoms->homenr,
1336 cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1337 if (useGpuPmeOnThisRank)
1339 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1340 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1343 // We need to copy coordinates when:
1344 // 1. Update is not offloaded
1345 // 2. The buffers were reinitialized on search step
1346 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1348 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1349 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1353 // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1354 // Otherwise the send will occur before the H2D coordinate transfer.
1355 if (!thisRankHasDuty(cr, DUTY_PME) && pmeSendCoordinatesFromGpu)
1357 /* Send particle coordinates to the pme nodes */
1358 gmx_pme_send_coordinates(fr,
1361 as_rvec_array(x.unpaddedArrayRef().data()),
1362 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1363 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1364 (stepWork.computeVirial || stepWork.computeEnergy),
1366 simulationWork.useGpuPmePpCommunication,
1367 reinitGpuPmePpComms,
1368 pmeSendCoordinatesFromGpu,
1369 localXReadyOnDevice,
1373 if (useGpuPmeOnThisRank)
1375 launchPmeGpuSpread(fr->pmedata,
1378 localXReadyOnDevice,
1379 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1383 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1385 /* do gridding for pair search */
1386 if (stepWork.doNeighborSearch)
1388 if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1390 fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1393 wallcycle_start(wcycle, ewcNS);
1394 if (!DOMAINDECOMP(cr))
1396 const rvec vzero = { 0.0_real, 0.0_real, 0.0_real };
1397 const rvec boxDiagonal = { box[XX][XX], box[YY][YY], box[ZZ][ZZ] };
1398 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1399 nbnxn_put_on_grid(nbv,
1405 { 0, mdatoms->homenr },
1408 x.unpaddedArrayRef(),
1411 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1415 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1416 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1417 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1420 nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1421 gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1424 wallcycle_stop(wcycle, ewcNS);
1426 /* initialize the GPU nbnxm atom data and bonded data structures */
1427 if (simulationWork.useGpuNonbonded)
1429 // Note: cycle counting only nononbondeds, gpuBonded counts internally
1430 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1431 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1432 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1433 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1434 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1438 /* Now we put all atoms on the grid, we can assign bonded
1439 * interactions to the GPU, where the grid order is
1440 * needed. Also the xq, f and fshift device buffers have
1441 * been reallocated if needed, so the bonded code can
1442 * learn about them. */
1443 // TODO the xq, f, and fshift buffers are now shared
1444 // resources, so they should be maintained by a
1445 // higher-level object than the nb module.
1446 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1448 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1449 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1450 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1454 // Need to run after the GPU-offload bonded interaction lists
1455 // are set up to be able to determine whether there is bonded work.
1456 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1457 inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1459 wallcycle_start_nocount(wcycle, ewcNS);
1460 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1461 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1462 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1464 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1466 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1467 wallcycle_stop(wcycle, ewcNS);
1469 if (stepWork.useGpuXBufferOps)
1471 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1474 if (simulationWork.useGpuBufferOps)
1476 setupGpuForceReductions(runScheduleWork, cr, fr);
1479 else if (!EI_TPI(inputrec.eI) && stepWork.computeNonbondedForces)
1481 if (stepWork.useGpuXBufferOps)
1483 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1484 nbv->convertCoordinatesGpu(AtomLocality::Local, stateGpu->getCoordinates(), localXReadyOnDevice);
1488 if (simulationWork.useGpuUpdate)
1490 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1491 GMX_ASSERT(haveCopiedXFromGpu,
1492 "a wait should only be triggered if copy has been scheduled");
1493 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1495 nbv->convertCoordinates(AtomLocality::Local, x.unpaddedArrayRef());
1499 if (simulationWork.useGpuNonbonded && (stepWork.computeNonbondedForces || domainWork.haveGpuBondedWork))
1501 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1503 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1504 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1505 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1506 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1508 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1510 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1511 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1512 // with X buffer ops offloaded to the GPU on all but the search steps
1514 // bonded work not split into separate local and non-local, so with DD
1515 // we can only launch the kernel after non-local coordinates have been received.
1516 if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1518 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1521 /* launch local nonbonded work on GPU */
1522 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1523 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1524 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1525 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1526 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1529 if (useGpuPmeOnThisRank)
1531 // In PME GPU and mixed mode we launch FFT / gather after the
1532 // X copy/transform to allow overlap as well as after the GPU NB
1533 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1534 // the nonbonded kernel.
1535 launchPmeGpuFftAndGather(fr->pmedata,
1536 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1541 /* Communicate coordinates and sum dipole if necessary +
1542 do non-local pair search */
1543 if (havePPDomainDecomposition(cr))
1545 if (stepWork.doNeighborSearch)
1547 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1548 wallcycle_start_nocount(wcycle, ewcNS);
1549 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1550 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1551 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1553 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1554 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1555 wallcycle_stop(wcycle, ewcNS);
1556 // TODO refactor this GPU halo exchange re-initialisation
1557 // to location in do_md where GPU halo exchange is
1558 // constructed at partitioning, after above stateGpu
1559 // re-initialization has similarly been refactored
1560 if (simulationWork.useGpuHaloExchange)
1562 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1567 if (stepWork.useGpuXHalo)
1569 // The following must be called after local setCoordinates (which records an event
1570 // when the coordinate data has been copied to the device).
1571 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1573 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1575 // non-local part of coordinate buffer must be copied back to host for CPU work
1576 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1581 if (simulationWork.useGpuUpdate)
1583 GMX_ASSERT(haveCopiedXFromGpu,
1584 "a wait should only be triggered if copy has been scheduled");
1585 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1587 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1590 if (stepWork.useGpuXBufferOps)
1592 if (!useGpuPmeOnThisRank && !stepWork.useGpuXHalo)
1594 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1596 nbv->convertCoordinatesGpu(AtomLocality::NonLocal,
1597 stateGpu->getCoordinates(),
1598 stateGpu->getCoordinatesReadyOnDeviceEvent(
1599 AtomLocality::NonLocal, simulationWork, stepWork));
1603 nbv->convertCoordinates(AtomLocality::NonLocal, x.unpaddedArrayRef());
1607 if (simulationWork.useGpuNonbonded)
1610 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1612 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1613 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1614 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1615 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1616 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1619 if (domainWork.haveGpuBondedWork)
1621 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1624 /* launch non-local nonbonded tasks on GPU */
1625 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1626 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1627 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1628 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1629 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1633 if (simulationWork.useGpuNonbonded && stepWork.computeNonbondedForces)
1635 /* launch D2H copy-back F */
1636 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1637 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1639 if (havePPDomainDecomposition(cr))
1641 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1643 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1644 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1646 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1648 fr->gpuBonded->launchEnergyTransfer();
1650 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1653 gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1654 if (fr->wholeMoleculeTransform)
1656 xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1659 DipoleData dipoleData;
1661 if (simulationWork.computeMuTot)
1663 const int start = 0;
1665 /* Calculate total (local) dipole moment in a temporary common array.
1666 * This makes it possible to sum them over nodes faster.
1668 gmx::ArrayRef<const gmx::RVec> xRef =
1669 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1673 gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1674 gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
1675 mdatoms->nChargePerturbed != 0,
1676 dipoleData.muStaging[0],
1677 dipoleData.muStaging[1]);
1679 reduceAndUpdateMuTot(
1680 &dipoleData, cr, (fr->efep != FreeEnergyPerturbationType::No), lambda, muTotal, ddBalanceRegionHandler);
1683 /* Reset energies */
1684 reset_enerdata(enerd);
1686 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1688 wallcycle_start(wcycle, ewcPPDURINGPME);
1689 dd_force_flop_start(cr->dd, nrnb);
1692 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1693 // this wait ensures that the D2H transfer is complete.
1694 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1695 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1697 GMX_ASSERT(haveCopiedXFromGpu, "a wait should only be triggered if copy has been scheduled");
1698 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1703 wallcycle_start(wcycle, ewcROT);
1704 do_rotation(cr, enforcedRotation, box, x.unpaddedConstArrayRef(), t, step, stepWork.doNeighborSearch);
1705 wallcycle_stop(wcycle, ewcROT);
1708 /* Start the force cycle counter.
1709 * Note that a different counter is used for dynamic load balancing.
1711 wallcycle_start(wcycle, ewcFORCE);
1713 /* Set up and clear force outputs:
1714 * forceOutMtsLevel0: everything except what is in the other two outputs
1715 * forceOutMtsLevel1: PME-mesh and listed-forces group 1
1716 * forceOutNonbonded: non-bonded forces
1717 * Without multiple time stepping all point to the same object.
1718 * With multiple time-stepping the use is different for MTS fast (level0 only) and slow steps.
1720 ForceOutputs forceOutMtsLevel0 = setupForceOutputs(
1721 &fr->forceHelperBuffers[0], force, domainWork, stepWork, havePPDomainDecomposition(cr), wcycle);
1723 // Force output for MTS combined forces, only set at level1 MTS steps
1724 std::optional<ForceOutputs> forceOutMts =
1725 (fr->useMts && stepWork.computeSlowForces)
1726 ? std::optional(setupForceOutputs(&fr->forceHelperBuffers[1],
1727 forceView->forceMtsCombinedWithPadding(),
1730 havePPDomainDecomposition(cr),
1734 ForceOutputs* forceOutMtsLevel1 =
1735 fr->useMts ? (stepWork.computeSlowForces ? &forceOutMts.value() : nullptr) : &forceOutMtsLevel0;
1737 const bool nonbondedAtMtsLevel1 = runScheduleWork->simulationWork.computeNonbondedAtMtsLevel1;
1739 ForceOutputs* forceOutNonbonded = nonbondedAtMtsLevel1 ? forceOutMtsLevel1 : &forceOutMtsLevel0;
1741 if (inputrec.bPull && pull_have_constraint(*pull_work))
1743 clear_pull_forces(pull_work);
1746 /* We calculate the non-bonded forces, when done on the CPU, here.
1747 * We do this before calling do_force_lowlevel, because in that
1748 * function, the listed forces are calculated before PME, which
1749 * does communication. With this order, non-bonded and listed
1750 * force calculation imbalance can be balanced out by the domain
1751 * decomposition load balancing.
1754 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1756 if (!useOrEmulateGpuNb)
1758 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1761 if (fr->efep != FreeEnergyPerturbationType::No && stepWork.computeNonbondedForces)
1763 /* Calculate the local and non-local free energy interactions here.
1764 * Happens here on the CPU both with and without GPU.
1766 nbv->dispatchFreeEnergyKernel(InteractionLocality::Local,
1768 x.unpaddedArrayRef(),
1769 &forceOutNonbonded->forceWithShiftForces(),
1770 gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1771 gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
1772 gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1773 gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr),
1774 inputrec.fepvals.get(),
1780 if (havePPDomainDecomposition(cr))
1782 nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal,
1784 x.unpaddedArrayRef(),
1785 &forceOutNonbonded->forceWithShiftForces(),
1786 gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1787 gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
1788 gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1789 gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr),
1790 inputrec.fepvals.get(),
1798 if (stepWork.computeNonbondedForces && !useOrEmulateGpuNb)
1800 if (havePPDomainDecomposition(cr))
1802 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1805 if (stepWork.computeForces)
1807 /* Add all the non-bonded force to the normal force array.
1808 * This can be split into a local and a non-local part when overlapping
1809 * communication with calculation with domain decomposition.
1811 wallcycle_stop(wcycle, ewcFORCE);
1812 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All,
1813 forceOutNonbonded->forceWithShiftForces().force());
1814 wallcycle_start_nocount(wcycle, ewcFORCE);
1817 /* If there are multiple fshift output buffers we need to reduce them */
1818 if (stepWork.computeVirial)
1820 /* This is not in a subcounter because it takes a
1821 negligible and constant-sized amount of time */
1822 nbnxn_atomdata_add_nbat_fshift_to_fshift(
1823 *nbv->nbat, forceOutNonbonded->forceWithShiftForces().shiftForces());
1827 // TODO Force flags should include haveFreeEnergyWork for this domain
1828 if (stepWork.useGpuXHalo && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1830 wallcycle_stop(wcycle, ewcFORCE);
1831 /* Wait for non-local coordinate data to be copied from device */
1832 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1833 wallcycle_start_nocount(wcycle, ewcFORCE);
1836 // Compute wall interactions, when present.
1837 // Note: should be moved to special forces.
1838 if (inputrec.nwall && stepWork.computeNonbondedForces)
1840 /* foreign lambda component for walls */
1841 real dvdl_walls = do_walls(inputrec,
1845 x.unpaddedConstArrayRef(),
1846 &forceOutMtsLevel0.forceWithVirial(),
1847 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1848 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
1850 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += dvdl_walls;
1853 if (stepWork.computeListedForces)
1855 /* Check whether we need to take into account PBC in listed interactions */
1856 bool needMolPbc = false;
1857 for (const auto& listedForces : fr->listedForces)
1859 if (listedForces.haveCpuListedForces(*fr->fcdata))
1861 needMolPbc = fr->bMolPBC;
1869 /* Since all atoms are in the rectangular or triclinic unit-cell,
1870 * only single box vector shifts (2 in x) are required.
1872 set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1875 for (int mtsIndex = 0; mtsIndex < (fr->useMts && stepWork.computeSlowForces ? 2 : 1); mtsIndex++)
1877 ListedForces& listedForces = fr->listedForces[mtsIndex];
1878 ForceOutputs& forceOut = (mtsIndex == 0 ? forceOutMtsLevel0 : *forceOutMtsLevel1);
1879 listedForces.calculate(wcycle,
1881 inputrec.fepvals.get(),
1895 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr,
1900 if (stepWork.computeSlowForces)
1902 calculateLongRangeNonbondeds(fr,
1908 x.unpaddedConstArrayRef(),
1909 &forceOutMtsLevel1->forceWithVirial(),
1913 dipoleData.muStateAB,
1915 ddBalanceRegionHandler);
1918 wallcycle_stop(wcycle, ewcFORCE);
1920 // VdW dispersion correction, only computed on master rank to avoid double counting
1921 if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1923 // Calculate long range corrections to pressure and energy
1924 const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
1925 box, lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)]);
1927 if (stepWork.computeEnergy)
1929 enerd->term[F_DISPCORR] = correction.energy;
1930 enerd->term[F_DVDL_VDW] += correction.dvdl;
1931 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += correction.dvdl;
1933 if (stepWork.computeVirial)
1935 correction.correctVirial(vir_force);
1936 enerd->term[F_PDISPCORR] = correction.pressure;
1940 computeSpecialForces(fplog,
1952 x.unpaddedArrayRef(),
1956 &forceOutMtsLevel0.forceWithVirial(),
1957 forceOutMtsLevel1 ? &forceOutMtsLevel1->forceWithVirial() : nullptr,
1960 stepWork.doNeighborSearch);
1962 if (havePPDomainDecomposition(cr) && stepWork.computeForces && stepWork.useGpuFHalo
1963 && domainWork.haveCpuLocalForceWork)
1965 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(), AtomLocality::Local);
1968 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
1969 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
1970 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFHalo),
1971 "The schedule below does not allow for nonbonded MTS with GPU halo exchange");
1972 // Will store the amount of cycles spent waiting for the GPU that
1973 // will be later used in the DLB accounting.
1974 float cycles_wait_gpu = 0;
1975 if (useOrEmulateGpuNb && stepWork.computeNonbondedForces)
1977 auto& forceWithShiftForces = forceOutNonbonded->forceWithShiftForces();
1979 /* wait for non-local forces (or calculate in emulation mode) */
1980 if (havePPDomainDecomposition(cr))
1982 if (simulationWork.useGpuNonbonded)
1984 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1987 AtomLocality::NonLocal,
1988 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
1989 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
1990 forceWithShiftForces.shiftForces(),
1995 wallcycle_start_nocount(wcycle, ewcFORCE);
1997 fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes, step, nrnb, wcycle);
1998 wallcycle_stop(wcycle, ewcFORCE);
2001 if (stepWork.useGpuFBufferOps)
2003 // TODO: move this into DomainLifetimeWorkload, including the second part of the
2004 // condition The bonded and free energy CPU tasks can have non-local force
2005 // contributions which are a dependency for the GPU force reduction.
2006 bool haveNonLocalForceContribInCpuBuffer =
2007 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
2009 if (haveNonLocalForceContribInCpuBuffer)
2011 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2012 AtomLocality::NonLocal);
2016 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->execute();
2018 if (!stepWork.useGpuFHalo)
2020 // copy from GPU input for dd_move_f()
2021 stateGpu->copyForcesFromGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2022 AtomLocality::NonLocal);
2027 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
2030 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
2032 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
2037 /* Combining the forces for multiple time stepping before the halo exchange, when possible,
2038 * avoids an extra halo exchange (when DD is used) and post-processing step.
2040 const bool combineMtsForcesBeforeHaloExchange =
2041 (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2042 && (legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0
2043 && !(stepWork.computeVirial || simulationWork.useGpuNonbonded || useGpuPmeOnThisRank));
2044 if (combineMtsForcesBeforeHaloExchange)
2046 const int numAtoms = havePPDomainDecomposition(cr) ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr;
2047 combineMtsForces(numAtoms,
2048 force.unpaddedArrayRef(),
2049 forceView->forceMtsCombined(),
2050 inputrec.mtsLevels[1].stepFactor);
2053 if (havePPDomainDecomposition(cr))
2055 /* We are done with the CPU compute.
2056 * We will now communicate the non-local forces.
2057 * If we use a GPU this will overlap with GPU work, so in that case
2058 * we do not close the DD force balancing region here.
2060 ddBalanceRegionHandler.closeAfterForceComputationCpu();
2062 if (stepWork.computeForces)
2065 if (stepWork.useGpuFHalo)
2067 // If there exist CPU forces, data from halo exchange should accumulate into these
2068 bool accumulateForces = domainWork.haveCpuLocalForceWork;
2069 if (!accumulateForces)
2071 // Force halo exchange will set a subset of local atoms with remote non-local data
2072 // First clear local portion of force array, so that untouched atoms are zero
2073 stateGpu->clearForcesOnGpu(AtomLocality::Local);
2075 communicateGpuHaloForces(*cr, accumulateForces);
2079 if (stepWork.useGpuFBufferOps)
2081 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
2084 // Without MTS or with MTS at slow steps with uncombined forces we need to
2085 // communicate the fast forces
2086 if (!fr->useMts || !combineMtsForcesBeforeHaloExchange)
2088 dd_move_f(cr->dd, &forceOutMtsLevel0.forceWithShiftForces(), wcycle);
2090 // With MTS we need to communicate the slow or combined (in forceOutMtsLevel1) forces
2091 if (fr->useMts && stepWork.computeSlowForces)
2093 dd_move_f(cr->dd, &forceOutMtsLevel1->forceWithShiftForces(), wcycle);
2099 // With both nonbonded and PME offloaded a GPU on the same rank, we use
2100 // an alternating wait/reduction scheme.
2101 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
2102 && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
2103 if (alternateGpuWait)
2105 alternatePmeNbGpuWaitReduce(fr->nbv.get(),
2110 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
2115 if (!alternateGpuWait && useGpuPmeOnThisRank)
2117 pme_gpu_wait_and_reduce(fr->pmedata,
2120 &forceOutMtsLevel1->forceWithVirial(),
2122 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]);
2125 /* Wait for local GPU NB outputs on the non-alternating wait path */
2126 if (!alternateGpuWait && stepWork.computeNonbondedForces && simulationWork.useGpuNonbonded)
2128 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
2129 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
2130 * but even with a step of 0.1 ms the difference is less than 1%
2133 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
2134 const float waitCycles = Nbnxm::gpu_wait_finish_task(
2137 AtomLocality::Local,
2138 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
2139 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
2140 forceOutNonbonded->forceWithShiftForces().shiftForces(),
2143 if (ddBalanceRegionHandler.useBalancingRegion())
2145 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
2146 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
2148 /* We measured few cycles, it could be that the kernel
2149 * and transfer finished earlier and there was no actual
2150 * wait time, only API call overhead.
2151 * Then the actual time could be anywhere between 0 and
2152 * cycles_wait_est. We will use half of cycles_wait_est.
2154 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
2156 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
2160 if (fr->nbv->emulateGpu())
2162 // NOTE: emulation kernel is not included in the balancing region,
2163 // but emulation mode does not target performance anyway
2164 wallcycle_start_nocount(wcycle, ewcFORCE);
2169 InteractionLocality::Local,
2170 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
2174 wallcycle_stop(wcycle, ewcFORCE);
2177 // If on GPU PME-PP comms path, receive forces from PME before GPU buffer ops
2178 // TODO refactor this and unify with below default-path call to the same function
2179 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces
2180 && simulationWork.useGpuPmePpCommunication)
2182 /* In case of node-splitting, the PP nodes receive the long-range
2183 * forces, virial and energy from the PME nodes here.
2185 pme_receive_force_ener(fr,
2187 &forceOutMtsLevel1->forceWithVirial(),
2189 simulationWork.useGpuPmePpCommunication,
2190 stepWork.useGpuPmeFReduction,
2195 /* Do the nonbonded GPU (or emulation) force buffer reduction
2196 * on the non-alternating path. */
2197 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
2198 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
2199 if (useOrEmulateGpuNb && !alternateGpuWait)
2201 if (stepWork.useGpuFBufferOps)
2203 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2205 // Flag to specify whether the CPU force buffer has contributions to
2206 // local atoms. This depends on whether there are CPU-based force tasks
2207 // or when DD is active the halo exchange has resulted in contributions
2208 // from the non-local part.
2209 const bool haveLocalForceContribInCpuBuffer =
2210 (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
2212 // TODO: move these steps as early as possible:
2213 // - CPU f H2D should be as soon as all CPU-side forces are done
2214 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
2215 // before the next CPU task that consumes the forces: vsite spread or update)
2216 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
2217 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
2218 // These should be unified.
2219 if (haveLocalForceContribInCpuBuffer && !stepWork.useGpuFHalo)
2221 // Note: AtomLocality::All is used for the non-DD case because, as in this
2222 // case copyForcesToGpu() uses a separate stream, it allows overlap of
2223 // CPU force H2D with GPU force tasks on all streams including those in the
2224 // local stream which would otherwise be implicit dependencies for the
2225 // transfer and would not overlap.
2226 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
2228 stateGpu->copyForcesToGpu(forceWithShift, locality);
2231 if (stepWork.computeNonbondedForces)
2233 fr->gpuForceReduction[gmx::AtomLocality::Local]->execute();
2236 // Copy forces to host if they are needed for update or if virtual sites are enabled.
2237 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
2238 // TODO: When the output flags will be included in step workload, this copy can be combined with the
2239 // copy call done in sim_utils(...) for the output.
2240 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
2241 // they should not be copied in do_md(...) for the output.
2242 if (!simulationWork.useGpuUpdate
2243 || (simulationWork.useGpuUpdate && DOMAINDECOMP(cr) && haveHostPmePpComms) || vsite)
2245 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
2246 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
2249 else if (stepWork.computeNonbondedForces)
2251 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2252 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
2256 launchGpuEndOfStepTasks(
2257 nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork, useGpuPmeOnThisRank, step, wcycle);
2259 if (DOMAINDECOMP(cr))
2261 dd_force_flop_stop(cr->dd, nrnb);
2264 const bool haveCombinedMtsForces = (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2265 && combineMtsForcesBeforeHaloExchange);
2266 if (stepWork.computeForces)
2268 postProcessForceWithShiftForces(
2269 nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutMtsLevel0, vir_force, *mdatoms, *fr, vsite, stepWork);
2271 if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2273 postProcessForceWithShiftForces(
2274 nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, *mdatoms, *fr, vsite, stepWork);
2278 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
2279 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
2280 && stepWork.computeSlowForces)
2282 /* In case of node-splitting, the PP nodes receive the long-range
2283 * forces, virial and energy from the PME nodes here.
2285 pme_receive_force_ener(fr,
2287 &forceOutMtsLevel1->forceWithVirial(),
2289 simulationWork.useGpuPmePpCommunication,
2294 if (stepWork.computeForces)
2296 /* If we don't use MTS or if we already combined the MTS forces before, we only
2297 * need to post-process one ForceOutputs object here, called forceOutCombined,
2298 * otherwise we have to post-process two outputs and then combine them.
2300 ForceOutputs& forceOutCombined = (haveCombinedMtsForces ? forceOutMts.value() : forceOutMtsLevel0);
2302 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutCombined, vir_force, mdatoms, fr, vsite, stepWork);
2304 if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2307 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, mdatoms, fr, vsite, stepWork);
2309 combineMtsForces(mdatoms->homenr,
2310 force.unpaddedArrayRef(),
2311 forceView->forceMtsCombined(),
2312 inputrec.mtsLevels[1].stepFactor);
2316 if (stepWork.computeEnergy)
2318 /* Compute the final potential energy terms */
2319 accumulatePotentialEnergies(enerd, lambda, inputrec.fepvals.get());
2321 if (!EI_TPI(inputrec.eI))
2323 checkPotentialEnergyValidity(step, *enerd, inputrec);
2327 /* In case we don't have constraints and are using GPUs, the next balancing
2328 * region starts here.
2329 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
2330 * virial calculation and COM pulling, is not thus not included in
2331 * the balance timing, which is ok as most tasks do communication.
2333 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);