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] +=
207 pull_potential(pull_work,
212 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Restraint)],
213 as_rvec_array(x.data()),
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(inputrec.pbcType,
696 rvec* f = as_rvec_array(forceWithVirialMtsLevel0->force_.data());
698 /* Add the forces from enforced rotation potentials (if any) */
701 wallcycle_start(wcycle, ewcROTadd);
702 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
703 wallcycle_stop(wcycle, ewcROTadd);
708 /* Note that since init_edsam() is called after the initialization
709 * of forcerec, edsam doesn't request the noVirSum force buffer.
710 * Thus if no other algorithm (e.g. PME) requires it, the forces
711 * here will contribute to the virial.
713 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
716 /* Add forces from interactive molecular dynamics (IMD), if any */
717 if (inputrec.bIMD && stepWork.computeForces)
719 imdSession->applyForces(f);
723 /*! \brief Launch the prepare_step and spread stages of PME GPU.
725 * \param[in] pmedata The PME structure
726 * \param[in] box The box matrix
727 * \param[in] stepWork Step schedule flags
728 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in the device memory.
729 * \param[in] lambdaQ The Coulomb lambda of the current state.
730 * \param[in] wcycle The wallcycle structure
732 static inline void launchPmeGpuSpread(gmx_pme_t* pmedata,
734 const StepWorkload& stepWork,
735 GpuEventSynchronizer* xReadyOnDevice,
737 gmx_wallcycle_t wcycle)
739 pme_gpu_prepare_computation(pmedata, box, wcycle, stepWork);
740 pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle, lambdaQ);
743 /*! \brief Launch the FFT and gather stages of PME GPU
745 * This function only implements setting the output forces (no accumulation).
747 * \param[in] pmedata The PME structure
748 * \param[in] lambdaQ The Coulomb lambda of the current system state.
749 * \param[in] wcycle The wallcycle structure
750 * \param[in] stepWork Step schedule flags
752 static void launchPmeGpuFftAndGather(gmx_pme_t* pmedata,
754 gmx_wallcycle_t wcycle,
755 const gmx::StepWorkload& stepWork)
757 pme_gpu_launch_complex_transforms(pmedata, wcycle, stepWork);
758 pme_gpu_launch_gather(pmedata, wcycle, lambdaQ);
762 * Polling wait for either of the PME or nonbonded GPU tasks.
764 * Instead of a static order in waiting for GPU tasks, this function
765 * polls checking which of the two tasks completes first, and does the
766 * associated force buffer reduction overlapped with the other task.
767 * By doing that, unlike static scheduling order, it can always overlap
768 * one of the reductions, regardless of the GPU task completion order.
770 * \param[in] nbv Nonbonded verlet structure
771 * \param[in,out] pmedata PME module data
772 * \param[in,out] forceOutputsNonbonded Force outputs for the non-bonded forces and shift forces
773 * \param[in,out] forceOutputsPme Force outputs for the PME forces and virial
774 * \param[in,out] enerd Energy data structure results are reduced into
775 * \param[in] lambdaQ The Coulomb lambda of the current system state.
776 * \param[in] stepWork Step schedule flags
777 * \param[in] wcycle The wallcycle structure
779 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
781 gmx::ForceOutputs* forceOutputsNonbonded,
782 gmx::ForceOutputs* forceOutputsPme,
783 gmx_enerdata_t* enerd,
785 const StepWorkload& stepWork,
786 gmx_wallcycle_t wcycle)
788 bool isPmeGpuDone = false;
789 bool isNbGpuDone = false;
791 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
793 while (!isPmeGpuDone || !isNbGpuDone)
797 GpuTaskCompletion completionType =
798 (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
799 isPmeGpuDone = pme_gpu_try_finish_task(
800 pmedata, stepWork, wcycle, &forceOutputsPme->forceWithVirial(), enerd, lambdaQ, completionType);
805 auto& forceBuffersNonbonded = forceOutputsNonbonded->forceWithShiftForces();
806 GpuTaskCompletion completionType =
807 (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
808 isNbGpuDone = Nbnxm::gpu_try_finish_task(
812 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
813 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
814 forceBuffersNonbonded.shiftForces(),
820 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceBuffersNonbonded.force());
826 /*! \brief Set up the different force buffers; also does clearing.
828 * \param[in] forceHelperBuffers Helper force buffers
829 * \param[in] force force array
830 * \param[in] stepWork Step schedule flags
831 * \param[out] wcycle wallcycle recording structure
833 * \returns Cleared force output structure
835 static ForceOutputs setupForceOutputs(ForceHelperBuffers* forceHelperBuffers,
836 gmx::ArrayRefWithPadding<gmx::RVec> force,
837 const StepWorkload& stepWork,
838 gmx_wallcycle_t wcycle)
840 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
842 /* NOTE: We assume fr->shiftForces is all zeros here */
843 gmx::ForceWithShiftForces forceWithShiftForces(
844 force, stepWork.computeVirial, forceHelperBuffers->shiftForces());
846 if (stepWork.computeForces)
848 /* Clear the short- and long-range forces */
849 clearRVecs(forceWithShiftForces.force(), true);
851 /* Clear the shift forces */
852 clearRVecs(forceWithShiftForces.shiftForces(), false);
855 /* If we need to compute the virial, we might need a separate
856 * force buffer for algorithms for which the virial is calculated
857 * directly, such as PME. Otherwise, forceWithVirial uses the
858 * the same force (f in legacy calls) buffer as other algorithms.
860 const bool useSeparateForceWithVirialBuffer =
861 (stepWork.computeForces
862 && (stepWork.computeVirial && forceHelperBuffers->haveDirectVirialContributions()));
863 /* forceWithVirial uses the local atom range only */
864 gmx::ForceWithVirial forceWithVirial(
865 useSeparateForceWithVirialBuffer ? forceHelperBuffers->forceBufferForDirectVirialContributions()
866 : force.unpaddedArrayRef(),
867 stepWork.computeVirial);
869 if (useSeparateForceWithVirialBuffer)
871 /* TODO: update comment
872 * We only compute forces on local atoms. Note that vsites can
873 * spread to non-local atoms, but that part of the buffer is
874 * cleared separately in the vsite spreading code.
876 clearRVecs(forceWithVirial.force_, true);
879 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
882 forceWithShiftForces, forceHelperBuffers->haveDirectVirialContributions(), forceWithVirial);
886 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
888 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec& inputrec,
889 const t_forcerec& fr,
890 const pull_t* pull_work,
892 const t_mdatoms& mdatoms,
893 const SimulationWorkload& simulationWork,
894 const StepWorkload& stepWork)
896 DomainLifetimeWorkload domainWork;
897 // Note that haveSpecialForces is constant over the whole run
898 domainWork.haveSpecialForces =
899 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
900 domainWork.haveCpuListedForceWork = false;
901 domainWork.haveCpuBondedWork = false;
902 for (const auto& listedForces : fr.listedForces)
904 if (listedForces.haveCpuListedForces(*fr.fcdata))
906 domainWork.haveCpuListedForceWork = true;
908 if (listedForces.haveCpuBondeds())
910 domainWork.haveCpuBondedWork = true;
913 domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
914 // Note that haveFreeEnergyWork is constant over the whole run
915 domainWork.haveFreeEnergyWork =
916 (fr.efep != FreeEnergyPerturbationType::No && mdatoms.nPerturbed != 0);
917 // We assume we have local force work if there are CPU
918 // force tasks including PME or nonbondeds.
919 domainWork.haveCpuLocalForceWork =
920 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
921 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
922 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
927 /*! \brief Set up force flag stuct from the force bitmask.
929 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
930 * \param[in] mtsLevels The multiple time-stepping levels, either empty or 2 levels
931 * \param[in] step The current MD step
932 * \param[in] simulationWork Simulation workload description.
933 * \param[in] rankHasPmeDuty If this rank computes PME.
935 * \returns New Stepworkload description.
937 static StepWorkload setupStepWorkload(const int legacyFlags,
938 ArrayRef<const gmx::MtsLevel> mtsLevels,
940 const SimulationWorkload& simulationWork,
941 const bool rankHasPmeDuty)
943 GMX_ASSERT(mtsLevels.empty() || mtsLevels.size() == 2, "Expect 0 or 2 MTS levels");
944 const bool computeSlowForces = (mtsLevels.empty() || step % mtsLevels[1].stepFactor == 0);
947 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
948 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
949 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
950 flags.computeSlowForces = computeSlowForces;
951 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
952 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
953 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
954 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
955 flags.computeNonbondedForces =
956 ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded
957 && !(simulationWork.computeNonbondedAtMtsLevel1 && !computeSlowForces);
958 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
960 if (simulationWork.useGpuBufferOps)
962 GMX_ASSERT(simulationWork.useGpuNonbonded,
963 "Can only offload buffer ops if nonbonded computation is also offloaded");
965 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
966 // on virial steps the CPU reduction path is taken
967 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
968 flags.useGpuPmeFReduction = flags.computeSlowForces && flags.useGpuFBufferOps && simulationWork.useGpuPme
969 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication);
970 flags.useGpuXHalo = simulationWork.useGpuHaloExchange;
971 flags.useGpuFHalo = simulationWork.useGpuHaloExchange && flags.useGpuFBufferOps;
977 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
979 * TODO: eliminate \p useGpuPmeOnThisRank when this is
980 * incorporated in DomainLifetimeWorkload.
982 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
983 gmx::GpuBonded* gpuBonded,
985 gmx_enerdata_t* enerd,
986 const gmx::MdrunScheduleWorkload& runScheduleWork,
987 bool useGpuPmeOnThisRank,
989 gmx_wallcycle_t wcycle)
991 if (runScheduleWork.simulationWork.useGpuNonbonded && runScheduleWork.stepWork.computeNonbondedForces)
993 /* Launch pruning before buffer clearing because the API overhead of the
994 * clear kernel launches can leave the GPU idle while it could be running
997 if (nbv->isDynamicPruningStepGpu(step))
999 nbv->dispatchPruneKernelGpu(step);
1002 /* now clear the GPU outputs while we finish the step on the CPU */
1003 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1004 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1005 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
1006 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1007 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1010 if (useGpuPmeOnThisRank)
1012 pme_gpu_reinit_computation(pmedata, wcycle);
1015 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
1017 // in principle this should be included in the DD balancing region,
1018 // but generally it is infrequent so we'll omit it for the sake of
1020 gpuBonded->waitAccumulateEnergyTerms(enerd);
1022 gpuBonded->clearEnergies();
1026 //! \brief Data structure to hold dipole-related data and staging arrays
1029 //! Dipole staging for fast summing over MPI
1030 gmx::DVec muStaging[2] = { { 0.0, 0.0, 0.0 } };
1031 //! Dipole staging for states A and B (index 0 and 1 resp.)
1032 gmx::RVec muStateAB[2] = { { 0.0_real, 0.0_real, 0.0_real } };
1036 static void reduceAndUpdateMuTot(DipoleData* dipoleData,
1037 const t_commrec* cr,
1038 const bool haveFreeEnergy,
1039 gmx::ArrayRef<const real> lambda,
1041 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1045 gmx_sumd(2 * DIM, dipoleData->muStaging[0], cr);
1046 ddBalanceRegionHandler.reopenRegionCpu();
1048 for (int i = 0; i < 2; i++)
1050 for (int j = 0; j < DIM; j++)
1052 dipoleData->muStateAB[i][j] = dipoleData->muStaging[i][j];
1056 if (!haveFreeEnergy)
1058 copy_rvec(dipoleData->muStateAB[0], muTotal);
1062 for (int j = 0; j < DIM; j++)
1064 muTotal[j] = (1.0 - lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)])
1065 * dipoleData->muStateAB[0][j]
1066 + lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]
1067 * dipoleData->muStateAB[1][j];
1072 /*! \brief Combines MTS level0 and level1 force buffes into a full and MTS-combined force buffer.
1074 * \param[in] numAtoms The number of atoms to combine forces for
1075 * \param[in,out] forceMtsLevel0 Input: F_level0, output: F_level0 + F_level1
1076 * \param[in,out] forceMts Input: F_level1, output: F_level0 + mtsFactor * F_level1
1077 * \param[in] mtsFactor The factor between the level0 and level1 time step
1079 static void combineMtsForces(const int numAtoms,
1080 ArrayRef<RVec> forceMtsLevel0,
1081 ArrayRef<RVec> forceMts,
1082 const real mtsFactor)
1084 const int gmx_unused numThreads = gmx_omp_nthreads_get(emntDefault);
1085 #pragma omp parallel for num_threads(numThreads) schedule(static)
1086 for (int i = 0; i < numAtoms; i++)
1088 const RVec forceMtsLevel0Tmp = forceMtsLevel0[i];
1089 forceMtsLevel0[i] += forceMts[i];
1090 forceMts[i] = forceMtsLevel0Tmp + mtsFactor * forceMts[i];
1094 /*! \brief Setup for the local and non-local GPU force reductions:
1095 * reinitialization plus the registration of forces and dependencies.
1097 * \param [in] runScheduleWork Schedule workload flag structure
1098 * \param [in] cr Communication record object
1099 * \param [in] fr Force record object
1101 static void setupGpuForceReductions(gmx::MdrunScheduleWorkload* runScheduleWork,
1102 const t_commrec* cr,
1106 nonbonded_verlet_t* nbv = fr->nbv.get();
1107 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1109 // (re-)initialize local GPU force reduction
1110 const bool accumulate =
1111 runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr);
1112 const int atomStart = 0;
1113 fr->gpuForceReduction[gmx::AtomLocality::Local]->reinit(stateGpu->getForces(),
1114 nbv->getNumAtoms(AtomLocality::Local),
1115 nbv->getGridIndices(),
1118 stateGpu->fReducedOnDevice());
1120 // register forces and add dependencies
1121 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerNbnxmForce(nbv->getGpuForces());
1123 if (runScheduleWork->simulationWork.useGpuPme
1124 && (thisRankHasDuty(cr, DUTY_PME) || runScheduleWork->simulationWork.useGpuPmePpCommunication))
1126 void* forcePtr = thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1127 : // PME force buffer on same GPU
1128 fr->pmePpCommGpu->getGpuForceStagingPtr(); // buffer received from other GPU
1129 fr->gpuForceReduction[gmx::AtomLocality::Local]->registerRvecForce(forcePtr);
1131 GpuEventSynchronizer* const pmeSynchronizer =
1132 (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1133 : // PME force buffer on same GPU
1134 fr->pmePpCommGpu->getForcesReadySynchronizer()); // buffer received from other GPU
1135 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(pmeSynchronizer);
1138 if ((runScheduleWork->domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr))
1139 && !runScheduleWork->simulationWork.useGpuHaloExchange)
1141 auto forcesReadyLocality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
1142 const bool useGpuForceBufferOps = true;
1143 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1144 stateGpu->getForcesReadyOnDeviceEvent(forcesReadyLocality, useGpuForceBufferOps));
1147 if (runScheduleWork->simulationWork.useGpuHaloExchange)
1149 fr->gpuForceReduction[gmx::AtomLocality::Local]->addDependency(
1150 cr->dd->gpuHaloExchange[0][0]->getForcesReadyOnDeviceEvent());
1153 if (havePPDomainDecomposition(cr))
1155 // (re-)initialize non-local GPU force reduction
1156 const bool accumulate = runScheduleWork->domainWork.haveCpuBondedWork
1157 || runScheduleWork->domainWork.haveFreeEnergyWork;
1158 const int atomStart = dd_numHomeAtoms(*cr->dd);
1159 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->reinit(stateGpu->getForces(),
1160 nbv->getNumAtoms(AtomLocality::NonLocal),
1161 nbv->getGridIndices(),
1165 // register forces and add dependencies
1166 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->registerNbnxmForce(nbv->getGpuForces());
1167 if (runScheduleWork->domainWork.haveCpuBondedWork || runScheduleWork->domainWork.haveFreeEnergyWork)
1169 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->addDependency(
1170 stateGpu->getForcesReadyOnDeviceEvent(AtomLocality::NonLocal, true));
1176 void do_force(FILE* fplog,
1177 const t_commrec* cr,
1178 const gmx_multisim_t* ms,
1179 const t_inputrec& inputrec,
1181 gmx_enfrot* enforcedRotation,
1182 gmx::ImdSession* imdSession,
1186 gmx_wallcycle_t wcycle,
1187 const gmx_localtop_t* top,
1189 gmx::ArrayRefWithPadding<gmx::RVec> x,
1190 const history_t* hist,
1191 gmx::ForceBuffersView* forceView,
1193 const t_mdatoms* mdatoms,
1194 gmx_enerdata_t* enerd,
1195 gmx::ArrayRef<const real> lambda,
1197 gmx::MdrunScheduleWorkload* runScheduleWork,
1198 gmx::VirtualSitesHandler* vsite,
1203 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1205 auto force = forceView->forceWithPadding();
1206 GMX_ASSERT(force.unpaddedArrayRef().ssize() >= fr->natoms_force_constr,
1207 "The size of the force buffer should be at least the number of atoms to compute "
1210 nonbonded_verlet_t* nbv = fr->nbv.get();
1211 interaction_const_t* ic = fr->ic.get();
1213 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
1215 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
1217 runScheduleWork->stepWork = setupStepWorkload(
1218 legacyFlags, inputrec.mtsLevels, step, simulationWork, thisRankHasDuty(cr, DUTY_PME));
1219 const StepWorkload& stepWork = runScheduleWork->stepWork;
1221 const bool useGpuPmeOnThisRank =
1222 simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces;
1224 /* At a search step we need to start the first balancing region
1225 * somewhere early inside the step after communication during domain
1226 * decomposition (and not during the previous step as usual).
1228 if (stepWork.doNeighborSearch)
1230 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
1233 clear_mat(vir_force);
1235 if (fr->pbcType != PbcType::No)
1237 /* Compute shift vectors every step,
1238 * because of pressure coupling or box deformation!
1240 if (stepWork.haveDynamicBox && stepWork.stateChanged)
1242 calc_shifts(box, fr->shift_vec);
1245 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
1246 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
1249 put_atoms_in_box_omp(fr->pbcType,
1251 x.unpaddedArrayRef().subArray(0, mdatoms->homenr),
1252 gmx_omp_nthreads_get(emntDefault));
1253 inc_nrnb(nrnb, eNR_SHIFTX, mdatoms->homenr);
1257 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1259 const bool pmeSendCoordinatesFromGpu =
1260 GMX_MPI && simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1261 const bool reinitGpuPmePpComms =
1262 GMX_MPI && simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1264 const auto localXReadyOnDevice = (useGpuPmeOnThisRank || simulationWork.useGpuBufferOps)
1265 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1266 AtomLocality::Local, simulationWork, stepWork)
1269 // Copy coordinate from the GPU if update is on the GPU and there
1270 // are forces to be computed on the CPU, or for the computation of
1271 // virial, or if host-side data will be transferred from this task
1272 // to a remote task for halo exchange or PME-PP communication. At
1273 // search steps the current coordinates are already on the host,
1274 // hence copy is not needed.
1275 const bool haveHostPmePpComms =
1276 !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1278 GMX_ASSERT(simulationWork.useGpuHaloExchange
1279 == ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange[0].empty())),
1280 "The GPU halo exchange is active, but it has not been constructed.");
1281 const bool haveHostHaloExchangeComms =
1282 havePPDomainDecomposition(cr) && !simulationWork.useGpuHaloExchange;
1284 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1285 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1286 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1287 || haveHostPmePpComms || haveHostHaloExchangeComms))
1289 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1290 haveCopiedXFromGpu = true;
1293 // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1294 // Otherwise the send will occur after H2D coordinate transfer.
1295 if (GMX_MPI && !thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu && stepWork.computeSlowForces)
1297 /* Send particle coordinates to the pme nodes */
1298 if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1300 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1303 gmx_pme_send_coordinates(fr,
1306 as_rvec_array(x.unpaddedArrayRef().data()),
1307 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1308 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1309 (stepWork.computeVirial || stepWork.computeEnergy),
1311 simulationWork.useGpuPmePpCommunication,
1312 reinitGpuPmePpComms,
1313 pmeSendCoordinatesFromGpu,
1314 localXReadyOnDevice,
1318 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1319 // The local coordinates can be copied right away.
1320 // NOTE: Consider moving this copy to right after they are updated and constrained,
1321 // if the later is not offloaded.
1322 if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1324 if (stepWork.doNeighborSearch)
1326 // TODO refactor this to do_md, after partitioning.
1327 stateGpu->reinit(mdatoms->homenr,
1328 cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1329 if (useGpuPmeOnThisRank)
1331 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1332 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1335 // We need to copy coordinates when:
1336 // 1. Update is not offloaded
1337 // 2. The buffers were reinitialized on search step
1338 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1340 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1341 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1345 // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1346 // Otherwise the send will occur before the H2D coordinate transfer.
1347 if (!thisRankHasDuty(cr, DUTY_PME) && pmeSendCoordinatesFromGpu)
1349 /* Send particle coordinates to the pme nodes */
1350 gmx_pme_send_coordinates(fr,
1353 as_rvec_array(x.unpaddedArrayRef().data()),
1354 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1355 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1356 (stepWork.computeVirial || stepWork.computeEnergy),
1358 simulationWork.useGpuPmePpCommunication,
1359 reinitGpuPmePpComms,
1360 pmeSendCoordinatesFromGpu,
1361 localXReadyOnDevice,
1365 if (useGpuPmeOnThisRank)
1367 launchPmeGpuSpread(fr->pmedata,
1370 localXReadyOnDevice,
1371 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1375 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1377 /* do gridding for pair search */
1378 if (stepWork.doNeighborSearch)
1380 if (fr->wholeMoleculeTransform && stepWork.stateChanged)
1382 fr->wholeMoleculeTransform->updateForAtomPbcJumps(x.unpaddedArrayRef(), box);
1385 wallcycle_start(wcycle, ewcNS);
1386 if (!DOMAINDECOMP(cr))
1388 const rvec vzero = { 0.0_real, 0.0_real, 0.0_real };
1389 const rvec boxDiagonal = { box[XX][XX], box[YY][YY], box[ZZ][ZZ] };
1390 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1391 nbnxn_put_on_grid(nbv,
1397 { 0, mdatoms->homenr },
1400 x.unpaddedArrayRef(),
1403 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1407 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1408 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1409 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1412 nbv->setAtomProperties(gmx::constArrayRefFromArray(mdatoms->typeA, mdatoms->nr),
1413 gmx::constArrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
1416 wallcycle_stop(wcycle, ewcNS);
1418 /* initialize the GPU nbnxm atom data and bonded data structures */
1419 if (simulationWork.useGpuNonbonded)
1421 // Note: cycle counting only nononbondeds, gpuBonded counts internally
1422 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1423 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1424 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1425 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1426 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1430 /* Now we put all atoms on the grid, we can assign bonded
1431 * interactions to the GPU, where the grid order is
1432 * needed. Also the xq, f and fshift device buffers have
1433 * been reallocated if needed, so the bonded code can
1434 * learn about them. */
1435 // TODO the xq, f, and fshift buffers are now shared
1436 // resources, so they should be maintained by a
1437 // higher-level object than the nb module.
1438 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1440 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1441 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1442 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1446 // Need to run after the GPU-offload bonded interaction lists
1447 // are set up to be able to determine whether there is bonded work.
1448 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1449 inputrec, *fr, pull_work, ed, *mdatoms, simulationWork, stepWork);
1451 wallcycle_start_nocount(wcycle, ewcNS);
1452 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1453 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1454 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1456 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1458 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1459 wallcycle_stop(wcycle, ewcNS);
1461 if (stepWork.useGpuXBufferOps)
1463 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1466 if (simulationWork.useGpuBufferOps)
1468 setupGpuForceReductions(runScheduleWork, cr, fr);
1471 else if (!EI_TPI(inputrec.eI) && stepWork.computeNonbondedForces)
1473 if (stepWork.useGpuXBufferOps)
1475 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1476 nbv->convertCoordinatesGpu(AtomLocality::Local, stateGpu->getCoordinates(), localXReadyOnDevice);
1480 if (simulationWork.useGpuUpdate)
1482 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1483 GMX_ASSERT(haveCopiedXFromGpu,
1484 "a wait should only be triggered if copy has been scheduled");
1485 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1487 nbv->convertCoordinates(AtomLocality::Local, x.unpaddedArrayRef());
1491 if (simulationWork.useGpuNonbonded && (stepWork.computeNonbondedForces || domainWork.haveGpuBondedWork))
1493 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1495 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1496 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1497 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1498 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1500 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1502 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1503 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1504 // with X buffer ops offloaded to the GPU on all but the search steps
1506 // bonded work not split into separate local and non-local, so with DD
1507 // we can only launch the kernel after non-local coordinates have been received.
1508 if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1510 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1513 /* launch local nonbonded work on GPU */
1514 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1515 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1516 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1517 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1518 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1521 if (useGpuPmeOnThisRank)
1523 // In PME GPU and mixed mode we launch FFT / gather after the
1524 // X copy/transform to allow overlap as well as after the GPU NB
1525 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1526 // the nonbonded kernel.
1527 launchPmeGpuFftAndGather(fr->pmedata,
1528 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
1533 /* Communicate coordinates and sum dipole if necessary +
1534 do non-local pair search */
1535 if (havePPDomainDecomposition(cr))
1537 if (stepWork.doNeighborSearch)
1539 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1540 wallcycle_start_nocount(wcycle, ewcNS);
1541 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1542 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1543 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1545 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1546 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1547 wallcycle_stop(wcycle, ewcNS);
1548 // TODO refactor this GPU halo exchange re-initialisation
1549 // to location in do_md where GPU halo exchange is
1550 // constructed at partitioning, after above stateGpu
1551 // re-initialization has similarly been refactored
1552 if (simulationWork.useGpuHaloExchange)
1554 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1559 if (stepWork.useGpuXHalo)
1561 // The following must be called after local setCoordinates (which records an event
1562 // when the coordinate data has been copied to the device).
1563 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1565 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1567 // non-local part of coordinate buffer must be copied back to host for CPU work
1568 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1573 if (simulationWork.useGpuUpdate)
1575 GMX_ASSERT(haveCopiedXFromGpu,
1576 "a wait should only be triggered if copy has been scheduled");
1577 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1579 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1582 if (stepWork.useGpuXBufferOps)
1584 if (!useGpuPmeOnThisRank && !stepWork.useGpuXHalo)
1586 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1588 nbv->convertCoordinatesGpu(AtomLocality::NonLocal,
1589 stateGpu->getCoordinates(),
1590 stateGpu->getCoordinatesReadyOnDeviceEvent(
1591 AtomLocality::NonLocal, simulationWork, stepWork));
1595 nbv->convertCoordinates(AtomLocality::NonLocal, x.unpaddedArrayRef());
1599 if (simulationWork.useGpuNonbonded)
1602 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1604 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1605 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1606 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1607 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1608 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1611 if (domainWork.haveGpuBondedWork)
1613 fr->gpuBonded->setPbcAndlaunchKernel(fr->pbcType, box, fr->bMolPBC, stepWork);
1616 /* launch non-local nonbonded tasks on GPU */
1617 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1618 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1619 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1620 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1621 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1625 if (simulationWork.useGpuNonbonded && stepWork.computeNonbondedForces)
1627 /* launch D2H copy-back F */
1628 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1629 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1631 if (havePPDomainDecomposition(cr))
1633 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1635 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1636 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1638 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1640 fr->gpuBonded->launchEnergyTransfer();
1642 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1645 gmx::ArrayRef<const gmx::RVec> xWholeMolecules;
1646 if (fr->wholeMoleculeTransform)
1648 xWholeMolecules = fr->wholeMoleculeTransform->wholeMoleculeCoordinates(x.unpaddedArrayRef(), box);
1651 DipoleData dipoleData;
1653 if (simulationWork.computeMuTot)
1655 const int start = 0;
1657 /* Calculate total (local) dipole moment in a temporary common array.
1658 * This makes it possible to sum them over nodes faster.
1660 gmx::ArrayRef<const gmx::RVec> xRef =
1661 (xWholeMolecules.empty() ? x.unpaddedArrayRef() : xWholeMolecules);
1667 mdatoms->nChargePerturbed,
1668 dipoleData.muStaging[0],
1669 dipoleData.muStaging[1]);
1671 reduceAndUpdateMuTot(
1672 &dipoleData, cr, (fr->efep != FreeEnergyPerturbationType::No), lambda, muTotal, ddBalanceRegionHandler);
1675 /* Reset energies */
1676 reset_enerdata(enerd);
1678 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1680 wallcycle_start(wcycle, ewcPPDURINGPME);
1681 dd_force_flop_start(cr->dd, nrnb);
1684 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1685 // this wait ensures that the D2H transfer is complete.
1686 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1687 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1689 GMX_ASSERT(haveCopiedXFromGpu, "a wait should only be triggered if copy has been scheduled");
1690 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1695 wallcycle_start(wcycle, ewcROT);
1696 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, stepWork.doNeighborSearch);
1697 wallcycle_stop(wcycle, ewcROT);
1700 /* Start the force cycle counter.
1701 * Note that a different counter is used for dynamic load balancing.
1703 wallcycle_start(wcycle, ewcFORCE);
1705 /* Set up and clear force outputs:
1706 * forceOutMtsLevel0: everything except what is in the other two outputs
1707 * forceOutMtsLevel1: PME-mesh and listed-forces group 1
1708 * forceOutNonbonded: non-bonded forces
1709 * Without multiple time stepping all point to the same object.
1710 * With multiple time-stepping the use is different for MTS fast (level0 only) and slow steps.
1712 ForceOutputs forceOutMtsLevel0 =
1713 setupForceOutputs(&fr->forceHelperBuffers[0], force, stepWork, wcycle);
1715 // Force output for MTS combined forces, only set at level1 MTS steps
1716 std::optional<ForceOutputs> forceOutMts =
1717 (fr->useMts && stepWork.computeSlowForces)
1718 ? std::optional(setupForceOutputs(&fr->forceHelperBuffers[1],
1719 forceView->forceMtsCombinedWithPadding(),
1724 ForceOutputs* forceOutMtsLevel1 =
1725 fr->useMts ? (stepWork.computeSlowForces ? &forceOutMts.value() : nullptr) : &forceOutMtsLevel0;
1727 const bool nonbondedAtMtsLevel1 = runScheduleWork->simulationWork.computeNonbondedAtMtsLevel1;
1729 ForceOutputs* forceOutNonbonded = nonbondedAtMtsLevel1 ? forceOutMtsLevel1 : &forceOutMtsLevel0;
1731 if (inputrec.bPull && pull_have_constraint(*pull_work))
1733 clear_pull_forces(pull_work);
1736 /* We calculate the non-bonded forces, when done on the CPU, here.
1737 * We do this before calling do_force_lowlevel, because in that
1738 * function, the listed forces are calculated before PME, which
1739 * does communication. With this order, non-bonded and listed
1740 * force calculation imbalance can be balanced out by the domain
1741 * decomposition load balancing.
1744 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1746 if (!useOrEmulateGpuNb)
1748 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1751 if (fr->efep != FreeEnergyPerturbationType::No && stepWork.computeNonbondedForces)
1753 /* Calculate the local and non-local free energy interactions here.
1754 * Happens here on the CPU both with and without GPU.
1756 nbv->dispatchFreeEnergyKernel(InteractionLocality::Local,
1758 as_rvec_array(x.unpaddedArrayRef().data()),
1759 &forceOutNonbonded->forceWithShiftForces(),
1761 inputrec.fepvals.get(),
1767 if (havePPDomainDecomposition(cr))
1769 nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal,
1771 as_rvec_array(x.unpaddedArrayRef().data()),
1772 &forceOutNonbonded->forceWithShiftForces(),
1774 inputrec.fepvals.get(),
1782 if (stepWork.computeNonbondedForces && !useOrEmulateGpuNb)
1784 if (havePPDomainDecomposition(cr))
1786 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step, nrnb, wcycle);
1789 if (stepWork.computeForces)
1791 /* Add all the non-bonded force to the normal force array.
1792 * This can be split into a local and a non-local part when overlapping
1793 * communication with calculation with domain decomposition.
1795 wallcycle_stop(wcycle, ewcFORCE);
1796 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All,
1797 forceOutNonbonded->forceWithShiftForces().force());
1798 wallcycle_start_nocount(wcycle, ewcFORCE);
1801 /* If there are multiple fshift output buffers we need to reduce them */
1802 if (stepWork.computeVirial)
1804 /* This is not in a subcounter because it takes a
1805 negligible and constant-sized amount of time */
1806 nbnxn_atomdata_add_nbat_fshift_to_fshift(
1807 *nbv->nbat, forceOutNonbonded->forceWithShiftForces().shiftForces());
1811 // TODO Force flags should include haveFreeEnergyWork for this domain
1812 if (stepWork.useGpuXHalo && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1814 wallcycle_stop(wcycle, ewcFORCE);
1815 /* Wait for non-local coordinate data to be copied from device */
1816 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1817 wallcycle_start_nocount(wcycle, ewcFORCE);
1820 // Compute wall interactions, when present.
1821 // Note: should be moved to special forces.
1822 if (inputrec.nwall && stepWork.computeNonbondedForces)
1824 /* foreign lambda component for walls */
1825 real dvdl_walls = do_walls(inputrec,
1829 x.unpaddedConstArrayRef(),
1830 &forceOutMtsLevel0.forceWithVirial(),
1831 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
1832 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
1834 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += dvdl_walls;
1837 if (stepWork.computeListedForces)
1839 /* Check whether we need to take into account PBC in listed interactions */
1840 bool needMolPbc = false;
1841 for (const auto& listedForces : fr->listedForces)
1843 if (listedForces.haveCpuListedForces(*fr->fcdata))
1845 needMolPbc = fr->bMolPBC;
1853 /* Since all atoms are in the rectangular or triclinic unit-cell,
1854 * only single box vector shifts (2 in x) are required.
1856 set_pbc_dd(&pbc, fr->pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
1859 for (int mtsIndex = 0; mtsIndex < (fr->useMts && stepWork.computeSlowForces ? 2 : 1); mtsIndex++)
1861 ListedForces& listedForces = fr->listedForces[mtsIndex];
1862 ForceOutputs& forceOut = (mtsIndex == 0 ? forceOutMtsLevel0 : *forceOutMtsLevel1);
1863 listedForces.calculate(wcycle,
1865 inputrec.fepvals.get(),
1879 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr,
1884 if (stepWork.computeSlowForces)
1886 calculateLongRangeNonbondeds(fr,
1892 x.unpaddedConstArrayRef(),
1893 &forceOutMtsLevel1->forceWithVirial(),
1897 dipoleData.muStateAB,
1899 ddBalanceRegionHandler);
1902 wallcycle_stop(wcycle, ewcFORCE);
1904 // VdW dispersion correction, only computed on master rank to avoid double counting
1905 if ((stepWork.computeEnergy || stepWork.computeVirial) && fr->dispersionCorrection && MASTER(cr))
1907 // Calculate long range corrections to pressure and energy
1908 const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
1909 box, lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)]);
1911 if (stepWork.computeEnergy)
1913 enerd->term[F_DISPCORR] = correction.energy;
1914 enerd->term[F_DVDL_VDW] += correction.dvdl;
1915 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] += correction.dvdl;
1917 if (stepWork.computeVirial)
1919 correction.correctVirial(vir_force);
1920 enerd->term[F_PDISPCORR] = correction.pressure;
1924 computeSpecialForces(fplog,
1936 x.unpaddedArrayRef(),
1940 &forceOutMtsLevel0.forceWithVirial(),
1941 forceOutMtsLevel1 ? &forceOutMtsLevel1->forceWithVirial() : nullptr,
1944 stepWork.doNeighborSearch);
1946 if (havePPDomainDecomposition(cr) && stepWork.computeForces && stepWork.useGpuFHalo
1947 && domainWork.haveCpuLocalForceWork)
1949 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(), AtomLocality::Local);
1952 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
1953 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
1954 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFHalo),
1955 "The schedule below does not allow for nonbonded MTS with GPU halo exchange");
1956 // Will store the amount of cycles spent waiting for the GPU that
1957 // will be later used in the DLB accounting.
1958 float cycles_wait_gpu = 0;
1959 if (useOrEmulateGpuNb && stepWork.computeNonbondedForces)
1961 auto& forceWithShiftForces = forceOutNonbonded->forceWithShiftForces();
1963 /* wait for non-local forces (or calculate in emulation mode) */
1964 if (havePPDomainDecomposition(cr))
1966 if (simulationWork.useGpuNonbonded)
1968 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1971 AtomLocality::NonLocal,
1972 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
1973 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
1974 forceWithShiftForces.shiftForces(),
1979 wallcycle_start_nocount(wcycle, ewcFORCE);
1981 fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes, step, nrnb, wcycle);
1982 wallcycle_stop(wcycle, ewcFORCE);
1985 if (stepWork.useGpuFBufferOps)
1987 // TODO: move this into DomainLifetimeWorkload, including the second part of the
1988 // condition The bonded and free energy CPU tasks can have non-local force
1989 // contributions which are a dependency for the GPU force reduction.
1990 bool haveNonLocalForceContribInCpuBuffer =
1991 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1993 if (haveNonLocalForceContribInCpuBuffer)
1995 stateGpu->copyForcesToGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
1996 AtomLocality::NonLocal);
2000 fr->gpuForceReduction[gmx::AtomLocality::NonLocal]->execute();
2002 if (!stepWork.useGpuFHalo)
2004 // copy from GPU input for dd_move_f()
2005 stateGpu->copyForcesFromGpu(forceOutMtsLevel0.forceWithShiftForces().force(),
2006 AtomLocality::NonLocal);
2011 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
2014 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
2016 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
2021 /* Combining the forces for multiple time stepping before the halo exchange, when possible,
2022 * avoids an extra halo exchange (when DD is used) and post-processing step.
2024 const bool combineMtsForcesBeforeHaloExchange =
2025 (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2026 && (legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0
2027 && !(stepWork.computeVirial || simulationWork.useGpuNonbonded || useGpuPmeOnThisRank));
2028 if (combineMtsForcesBeforeHaloExchange)
2030 const int numAtoms = havePPDomainDecomposition(cr) ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr;
2031 combineMtsForces(numAtoms,
2032 force.unpaddedArrayRef(),
2033 forceView->forceMtsCombined(),
2034 inputrec.mtsLevels[1].stepFactor);
2037 if (havePPDomainDecomposition(cr))
2039 /* We are done with the CPU compute.
2040 * We will now communicate the non-local forces.
2041 * If we use a GPU this will overlap with GPU work, so in that case
2042 * we do not close the DD force balancing region here.
2044 ddBalanceRegionHandler.closeAfterForceComputationCpu();
2046 if (stepWork.computeForces)
2049 if (stepWork.useGpuFHalo)
2051 // If there exist CPU forces, data from halo exchange should accumulate into these
2052 bool accumulateForces = domainWork.haveCpuLocalForceWork;
2053 if (!accumulateForces)
2055 // Force halo exchange will set a subset of local atoms with remote non-local data
2056 // First clear local portion of force array, so that untouched atoms are zero
2057 stateGpu->clearForcesOnGpu(AtomLocality::Local);
2059 communicateGpuHaloForces(*cr, accumulateForces);
2063 if (stepWork.useGpuFBufferOps)
2065 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
2068 // Without MTS or with MTS at slow steps with uncombined forces we need to
2069 // communicate the fast forces
2070 if (!fr->useMts || !combineMtsForcesBeforeHaloExchange)
2072 dd_move_f(cr->dd, &forceOutMtsLevel0.forceWithShiftForces(), wcycle);
2074 // With MTS we need to communicate the slow or combined (in forceOutMtsLevel1) forces
2075 if (fr->useMts && stepWork.computeSlowForces)
2077 dd_move_f(cr->dd, &forceOutMtsLevel1->forceWithShiftForces(), wcycle);
2083 // With both nonbonded and PME offloaded a GPU on the same rank, we use
2084 // an alternating wait/reduction scheme.
2085 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
2086 && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
2087 if (alternateGpuWait)
2089 alternatePmeNbGpuWaitReduce(fr->nbv.get(),
2094 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
2099 if (!alternateGpuWait && useGpuPmeOnThisRank)
2101 pme_gpu_wait_and_reduce(fr->pmedata,
2104 &forceOutMtsLevel1->forceWithVirial(),
2106 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)]);
2109 /* Wait for local GPU NB outputs on the non-alternating wait path */
2110 if (!alternateGpuWait && stepWork.computeNonbondedForces && simulationWork.useGpuNonbonded)
2112 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
2113 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
2114 * but even with a step of 0.1 ms the difference is less than 1%
2117 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
2118 const float waitCycles = Nbnxm::gpu_wait_finish_task(
2121 AtomLocality::Local,
2122 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
2123 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
2124 forceOutNonbonded->forceWithShiftForces().shiftForces(),
2127 if (ddBalanceRegionHandler.useBalancingRegion())
2129 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
2130 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
2132 /* We measured few cycles, it could be that the kernel
2133 * and transfer finished earlier and there was no actual
2134 * wait time, only API call overhead.
2135 * Then the actual time could be anywhere between 0 and
2136 * cycles_wait_est. We will use half of cycles_wait_est.
2138 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
2140 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
2144 if (fr->nbv->emulateGpu())
2146 // NOTE: emulation kernel is not included in the balancing region,
2147 // but emulation mode does not target performance anyway
2148 wallcycle_start_nocount(wcycle, ewcFORCE);
2153 InteractionLocality::Local,
2154 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
2158 wallcycle_stop(wcycle, ewcFORCE);
2161 // If on GPU PME-PP comms path, receive forces from PME before GPU buffer ops
2162 // TODO refactor this and unify with below default-path call to the same function
2163 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && stepWork.computeSlowForces
2164 && simulationWork.useGpuPmePpCommunication)
2166 /* In case of node-splitting, the PP nodes receive the long-range
2167 * forces, virial and energy from the PME nodes here.
2169 pme_receive_force_ener(fr,
2171 &forceOutMtsLevel1->forceWithVirial(),
2173 simulationWork.useGpuPmePpCommunication,
2174 stepWork.useGpuPmeFReduction,
2179 /* Do the nonbonded GPU (or emulation) force buffer reduction
2180 * on the non-alternating path. */
2181 GMX_ASSERT(!(nonbondedAtMtsLevel1 && stepWork.useGpuFBufferOps),
2182 "The schedule below does not allow for nonbonded MTS with GPU buffer ops");
2183 if (useOrEmulateGpuNb && !alternateGpuWait)
2185 if (stepWork.useGpuFBufferOps)
2187 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2189 // Flag to specify whether the CPU force buffer has contributions to
2190 // local atoms. This depends on whether there are CPU-based force tasks
2191 // or when DD is active the halo exchange has resulted in contributions
2192 // from the non-local part.
2193 const bool haveLocalForceContribInCpuBuffer =
2194 (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
2196 // TODO: move these steps as early as possible:
2197 // - CPU f H2D should be as soon as all CPU-side forces are done
2198 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
2199 // before the next CPU task that consumes the forces: vsite spread or update)
2200 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
2201 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
2202 // These should be unified.
2203 if (haveLocalForceContribInCpuBuffer && !stepWork.useGpuFHalo)
2205 // Note: AtomLocality::All is used for the non-DD case because, as in this
2206 // case copyForcesToGpu() uses a separate stream, it allows overlap of
2207 // CPU force H2D with GPU force tasks on all streams including those in the
2208 // local stream which would otherwise be implicit dependencies for the
2209 // transfer and would not overlap.
2210 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
2212 stateGpu->copyForcesToGpu(forceWithShift, locality);
2215 if (stepWork.computeNonbondedForces)
2217 fr->gpuForceReduction[gmx::AtomLocality::Local]->execute();
2220 // Copy forces to host if they are needed for update or if virtual sites are enabled.
2221 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
2222 // TODO: When the output flags will be included in step workload, this copy can be combined with the
2223 // copy call done in sim_utils(...) for the output.
2224 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
2225 // they should not be copied in do_md(...) for the output.
2226 if (!simulationWork.useGpuUpdate
2227 || (simulationWork.useGpuUpdate && DOMAINDECOMP(cr) && haveHostPmePpComms) || vsite)
2229 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
2230 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
2233 else if (stepWork.computeNonbondedForces)
2235 ArrayRef<gmx::RVec> forceWithShift = forceOutNonbonded->forceWithShiftForces().force();
2236 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
2240 launchGpuEndOfStepTasks(
2241 nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork, useGpuPmeOnThisRank, step, wcycle);
2243 if (DOMAINDECOMP(cr))
2245 dd_force_flop_stop(cr->dd, nrnb);
2248 const bool haveCombinedMtsForces = (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
2249 && combineMtsForcesBeforeHaloExchange);
2250 if (stepWork.computeForces)
2252 postProcessForceWithShiftForces(
2253 nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutMtsLevel0, vir_force, *mdatoms, *fr, vsite, stepWork);
2255 if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2257 postProcessForceWithShiftForces(
2258 nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, *mdatoms, *fr, vsite, stepWork);
2262 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
2263 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
2264 && stepWork.computeSlowForces)
2266 /* In case of node-splitting, the PP nodes receive the long-range
2267 * forces, virial and energy from the PME nodes here.
2269 pme_receive_force_ener(fr,
2271 &forceOutMtsLevel1->forceWithVirial(),
2273 simulationWork.useGpuPmePpCommunication,
2278 if (stepWork.computeForces)
2280 /* If we don't use MTS or if we already combined the MTS forces before, we only
2281 * need to post-process one ForceOutputs object here, called forceOutCombined,
2282 * otherwise we have to post-process two outputs and then combine them.
2284 ForceOutputs& forceOutCombined = (haveCombinedMtsForces ? forceOutMts.value() : forceOutMtsLevel0);
2286 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), &forceOutCombined, vir_force, mdatoms, fr, vsite, stepWork);
2288 if (fr->useMts && stepWork.computeSlowForces && !haveCombinedMtsForces)
2291 cr, step, nrnb, wcycle, box, x.unpaddedArrayRef(), forceOutMtsLevel1, vir_force, mdatoms, fr, vsite, stepWork);
2293 combineMtsForces(mdatoms->homenr,
2294 force.unpaddedArrayRef(),
2295 forceView->forceMtsCombined(),
2296 inputrec.mtsLevels[1].stepFactor);
2300 if (stepWork.computeEnergy)
2302 /* Compute the final potential energy terms */
2303 accumulatePotentialEnergies(enerd, lambda, inputrec.fepvals.get());
2305 if (!EI_TPI(inputrec.eI))
2307 checkPotentialEnergyValidity(step, *enerd, inputrec);
2311 /* In case we don't have constraints and are using GPUs, the next balancing
2312 * region starts here.
2313 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
2314 * virial calculation and COM pulling, is not thus not included in
2315 * the balance timing, which is ok as most tasks do communication.
2317 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);