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,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "gromacs/awh/awh.h"
49 #include "gromacs/domdec/dlbtiming.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/domdec/partition.h"
53 #include "gromacs/essentialdynamics/edsam.h"
54 #include "gromacs/ewald/pme.h"
55 #include "gromacs/gmxlib/chargegroup.h"
56 #include "gromacs/gmxlib/network.h"
57 #include "gromacs/gmxlib/nrnb.h"
58 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
59 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
60 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
61 #include "gromacs/gpu_utils/gpu_utils.h"
62 #include "gromacs/imd/imd.h"
63 #include "gromacs/listed_forces/bonded.h"
64 #include "gromacs/listed_forces/disre.h"
65 #include "gromacs/listed_forces/gpubonded.h"
66 #include "gromacs/listed_forces/listed_forces.h"
67 #include "gromacs/listed_forces/manage_threading.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/enerdata_utils.h"
78 #include "gromacs/mdlib/force.h"
79 #include "gromacs/mdlib/forcerec.h"
80 #include "gromacs/mdlib/gmx_omp_nthreads.h"
81 #include "gromacs/mdlib/ppforceworkload.h"
82 #include "gromacs/mdlib/qmmm.h"
83 #include "gromacs/mdlib/update.h"
84 #include "gromacs/mdtypes/commrec.h"
85 #include "gromacs/mdtypes/enerdata.h"
86 #include "gromacs/mdtypes/forceoutput.h"
87 #include "gromacs/mdtypes/iforceprovider.h"
88 #include "gromacs/mdtypes/inputrec.h"
89 #include "gromacs/mdtypes/md_enums.h"
90 #include "gromacs/mdtypes/state.h"
91 #include "gromacs/nbnxm/atomdata.h"
92 #include "gromacs/nbnxm/gpu_data_mgmt.h"
93 #include "gromacs/nbnxm/nbnxm.h"
94 #include "gromacs/pbcutil/ishift.h"
95 #include "gromacs/pbcutil/mshift.h"
96 #include "gromacs/pbcutil/pbc.h"
97 #include "gromacs/pulling/pull.h"
98 #include "gromacs/pulling/pull_rotation.h"
99 #include "gromacs/timing/cyclecounter.h"
100 #include "gromacs/timing/gpu_timing.h"
101 #include "gromacs/timing/wallcycle.h"
102 #include "gromacs/timing/wallcyclereporting.h"
103 #include "gromacs/timing/walltime_accounting.h"
104 #include "gromacs/topology/topology.h"
105 #include "gromacs/utility/arrayref.h"
106 #include "gromacs/utility/basedefinitions.h"
107 #include "gromacs/utility/cstringutil.h"
108 #include "gromacs/utility/exceptions.h"
109 #include "gromacs/utility/fatalerror.h"
110 #include "gromacs/utility/gmxassert.h"
111 #include "gromacs/utility/gmxmpi.h"
112 #include "gromacs/utility/logger.h"
113 #include "gromacs/utility/smalloc.h"
114 #include "gromacs/utility/strconvert.h"
115 #include "gromacs/utility/sysinfo.h"
117 // TODO: this environment variable allows us to verify before release
118 // that on less common architectures the total cost of polling is not larger than
119 // a blocking wait (so polling does not introduce overhead when the static
120 // PME-first ordering would suffice).
121 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
123 // environment variable to enable GPU buffer ops, to alow incremental and optional
124 // introduction of this functionality.
125 // TODO eventially tie this in with other existing GPU flags.
126 static const bool c_enableGpuBufOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
130 static void sum_forces(rvec f[], gmx::ArrayRef<const gmx::RVec> forceToAdd)
132 const int end = forceToAdd.size();
134 int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
135 #pragma omp parallel for num_threads(nt) schedule(static)
136 for (int i = 0; i < end; i++)
138 rvec_inc(f[i], forceToAdd[i]);
142 static void calc_virial(int start, int homenr, const rvec x[], const rvec f[],
143 tensor vir_part, const t_graph *graph, const matrix box,
144 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
146 /* The short-range virial from surrounding boxes */
147 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
148 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
150 /* Calculate partial virial, for local atoms only, based on short range.
151 * Total virial is computed in global_stat, called from do_md
153 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
154 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
158 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
162 static void pull_potential_wrapper(const t_commrec *cr,
163 const t_inputrec *ir,
164 const matrix box, gmx::ArrayRef<const gmx::RVec> x,
165 gmx::ForceWithVirial *force,
166 const t_mdatoms *mdatoms,
167 gmx_enerdata_t *enerd,
171 gmx_wallcycle_t wcycle)
176 /* Calculate the center of mass forces, this requires communication,
177 * which is why pull_potential is called close to other communication.
179 wallcycle_start(wcycle, ewcPULLPOT);
180 set_pbc(&pbc, ir->ePBC, box);
182 enerd->term[F_COM_PULL] +=
183 pull_potential(pull_work, mdatoms, &pbc,
184 cr, t, lambda[efptRESTRAINT], as_rvec_array(x.data()), force, &dvdl);
185 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
186 wallcycle_stop(wcycle, ewcPULLPOT);
189 static void pme_receive_force_ener(const t_commrec *cr,
190 gmx::ForceWithVirial *forceWithVirial,
191 gmx_enerdata_t *enerd,
192 gmx_wallcycle_t wcycle)
194 real e_q, e_lj, dvdl_q, dvdl_lj;
195 float cycles_ppdpme, cycles_seppme;
197 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
198 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
200 /* In case of node-splitting, the PP nodes receive the long-range
201 * forces, virial and energy from the PME nodes here.
203 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
206 gmx_pme_receive_f(cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
208 enerd->term[F_COUL_RECIP] += e_q;
209 enerd->term[F_LJ_RECIP] += e_lj;
210 enerd->dvdl_lin[efptCOUL] += dvdl_q;
211 enerd->dvdl_lin[efptVDW] += dvdl_lj;
215 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
217 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
220 static void print_large_forces(FILE *fp,
228 real force2Tolerance = gmx::square(forceTolerance);
229 gmx::index numNonFinite = 0;
230 for (int i = 0; i < md->homenr; i++)
232 real force2 = norm2(f[i]);
233 bool nonFinite = !std::isfinite(force2);
234 if (force2 >= force2Tolerance || nonFinite)
236 fprintf(fp, "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
238 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
245 if (numNonFinite > 0)
247 /* Note that with MPI this fatal call on one rank might interrupt
248 * the printing on other ranks. But we can only avoid that with
249 * an expensive MPI barrier that we would need at each step.
251 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
255 static void post_process_forces(const t_commrec *cr,
258 gmx_wallcycle_t wcycle,
259 const gmx_localtop_t *top,
263 gmx::ForceWithVirial *forceWithVirial,
265 const t_mdatoms *mdatoms,
266 const t_graph *graph,
267 const t_forcerec *fr,
268 const gmx_vsite_t *vsite,
271 if (fr->haveDirectVirialContributions)
273 rvec *fDirectVir = as_rvec_array(forceWithVirial->force_.data());
277 /* Spread the mesh force on virtual sites to the other particles...
278 * This is parallellized. MPI communication is performed
279 * if the constructing atoms aren't local.
281 matrix virial = { { 0 } };
282 spread_vsite_f(vsite, x, fDirectVir, nullptr,
283 (flags & GMX_FORCE_VIRIAL) != 0, virial,
285 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
286 forceWithVirial->addVirialContribution(virial);
289 if (flags & GMX_FORCE_VIRIAL)
291 /* Now add the forces, this is local */
292 sum_forces(f, forceWithVirial->force_);
294 /* Add the direct virial contributions */
295 GMX_ASSERT(forceWithVirial->computeVirial_, "forceWithVirial should request virial computation when we request the virial");
296 m_add(vir_force, forceWithVirial->getVirial(), vir_force);
300 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
305 if (fr->print_force >= 0)
307 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
311 static void do_nb_verlet(t_forcerec *fr,
312 const interaction_const_t *ic,
313 gmx_enerdata_t *enerd,
315 const Nbnxm::InteractionLocality ilocality,
319 gmx_wallcycle_t wcycle)
321 if (!(flags & GMX_FORCE_NONBONDED))
323 /* skip non-bonded calculation */
327 nonbonded_verlet_t *nbv = fr->nbv.get();
329 /* GPU kernel launch overhead is already timed separately */
330 if (fr->cutoff_scheme != ecutsVERLET)
332 gmx_incons("Invalid cut-off scheme passed!");
337 /* When dynamic pair-list pruning is requested, we need to prune
338 * at nstlistPrune steps.
340 if (nbv->isDynamicPruningStepCpu(step))
342 /* Prune the pair-list beyond fr->ic->rlistPrune using
343 * the current coordinates of the atoms.
345 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
346 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
347 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
350 wallcycle_sub_start(wcycle, ewcsNONBONDED);
353 nbv->dispatchNonbondedKernel(ilocality, *ic, flags, clearF, *fr, enerd, nrnb, wcycle);
357 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
361 static inline void clear_rvecs_omp(int n, rvec v[])
363 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, n);
365 /* Note that we would like to avoid this conditional by putting it
366 * into the omp pragma instead, but then we still take the full
367 * omp parallel for overhead (at least with gcc5).
371 for (int i = 0; i < n; i++)
378 #pragma omp parallel for num_threads(nth) schedule(static)
379 for (int i = 0; i < n; i++)
386 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
388 * \param groupOptions Group options, containing T-coupling options
390 static real averageKineticEnergyEstimate(const t_grpopts &groupOptions)
392 real nrdfCoupled = 0;
393 real nrdfUncoupled = 0;
394 real kineticEnergy = 0;
395 for (int g = 0; g < groupOptions.ngtc; g++)
397 if (groupOptions.tau_t[g] >= 0)
399 nrdfCoupled += groupOptions.nrdf[g];
400 kineticEnergy += groupOptions.nrdf[g]*0.5*groupOptions.ref_t[g]*BOLTZ;
404 nrdfUncoupled += groupOptions.nrdf[g];
408 /* This conditional with > also catches nrdf=0 */
409 if (nrdfCoupled > nrdfUncoupled)
411 return kineticEnergy*(nrdfCoupled + nrdfUncoupled)/nrdfCoupled;
419 /*! \brief This routine checks that the potential energy is finite.
421 * Always checks that the potential energy is finite. If step equals
422 * inputrec.init_step also checks that the magnitude of the potential energy
423 * is reasonable. Terminates with a fatal error when a check fails.
424 * Note that passing this check does not guarantee finite forces,
425 * since those use slightly different arithmetics. But in most cases
426 * there is just a narrow coordinate range where forces are not finite
427 * and energies are finite.
429 * \param[in] step The step number, used for checking and printing
430 * \param[in] enerd The energy data; the non-bonded group energies need to be added to enerd.term[F_EPOT] before calling this routine
431 * \param[in] inputrec The input record
433 static void checkPotentialEnergyValidity(int64_t step,
434 const gmx_enerdata_t &enerd,
435 const t_inputrec &inputrec)
437 /* Threshold valid for comparing absolute potential energy against
438 * the kinetic energy. Normally one should not consider absolute
439 * potential energy values, but with a factor of one million
440 * we should never get false positives.
442 constexpr real c_thresholdFactor = 1e6;
444 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
445 real averageKineticEnergy = 0;
446 /* We only check for large potential energy at the initial step,
447 * because that is by far the most likely step for this too occur
448 * and because computing the average kinetic energy is not free.
449 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
450 * before they become NaN.
452 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
454 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
457 if (energyIsNotFinite || (averageKineticEnergy > 0 &&
458 enerd.term[F_EPOT] > c_thresholdFactor*averageKineticEnergy))
460 gmx_fatal(FARGS, "Step %" PRId64 ": The total potential energy is %g, which is %s. The LJ and electrostatic contributions to the energy are %g and %g, respectively. A %s potential energy can be caused by overlapping interactions in bonded interactions or very large%s coordinate values. Usually this is caused by a badly- or non-equilibrated initial configuration, incorrect interactions or parameters in the topology.",
463 energyIsNotFinite ? "not finite" : "extremely high",
465 enerd.term[F_COUL_SR],
466 energyIsNotFinite ? "non-finite" : "very high",
467 energyIsNotFinite ? " or Nan" : "");
471 /*! \brief Return true if there are special forces computed this step.
473 * The conditionals exactly correspond to those in computeSpecialForces().
476 haveSpecialForces(const t_inputrec *inputrec,
477 ForceProviders *forceProviders,
478 const pull_t *pull_work,
482 const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
485 ((computeForces && forceProviders->hasForceProvider()) || // forceProviders
486 (inputrec->bPull && pull_have_potential(pull_work)) || // pull
487 inputrec->bRot || // enforced rotation
488 (ed != nullptr) || // flooding
489 (inputrec->bIMD && computeForces)); // IMD
492 /*! \brief Compute forces and/or energies for special algorithms
494 * The intention is to collect all calls to algorithms that compute
495 * forces on local atoms only and that do not contribute to the local
496 * virial sum (but add their virial contribution separately).
497 * Eventually these should likely all become ForceProviders.
498 * Within this function the intention is to have algorithms that do
499 * global communication at the end, so global barriers within the MD loop
500 * are as close together as possible.
502 * \param[in] fplog The log file
503 * \param[in] cr The communication record
504 * \param[in] inputrec The input record
505 * \param[in] awh The Awh module (nullptr if none in use).
506 * \param[in] enforcedRotation Enforced rotation module.
507 * \param[in] imdSession The IMD session
508 * \param[in] pull_work The pull work structure.
509 * \param[in] step The current MD step
510 * \param[in] t The current time
511 * \param[in,out] wcycle Wallcycle accounting struct
512 * \param[in,out] forceProviders Pointer to a list of force providers
513 * \param[in] box The unit cell
514 * \param[in] x The coordinates
515 * \param[in] mdatoms Per atom properties
516 * \param[in] lambda Array of free-energy lambda values
517 * \param[in] forceFlags Flags that tell whether we should compute forces/energies/virial
518 * \param[in,out] forceWithVirial Force and virial buffers
519 * \param[in,out] enerd Energy buffer
520 * \param[in,out] ed Essential dynamics pointer
521 * \param[in] bNS Tells if we did neighbor searching this step, used for ED sampling
523 * \todo Remove bNS, which is used incorrectly.
524 * \todo Convert all other algorithms called here to ForceProviders.
527 computeSpecialForces(FILE *fplog,
529 const t_inputrec *inputrec,
531 gmx_enfrot *enforcedRotation,
532 gmx::ImdSession *imdSession,
536 gmx_wallcycle_t wcycle,
537 ForceProviders *forceProviders,
539 gmx::ArrayRef<const gmx::RVec> x,
540 const t_mdatoms *mdatoms,
543 gmx::ForceWithVirial *forceWithVirial,
544 gmx_enerdata_t *enerd,
548 const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
550 /* NOTE: Currently all ForceProviders only provide forces.
551 * When they also provide energies, remove this conditional.
555 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
556 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
558 /* Collect forces from modules */
559 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
562 if (inputrec->bPull && pull_have_potential(pull_work))
564 pull_potential_wrapper(cr, inputrec, box, x,
566 mdatoms, enerd, pull_work, lambda, t,
571 enerd->term[F_COM_PULL] +=
572 awh->applyBiasForcesAndUpdateBias(inputrec->ePBC, *mdatoms, box,
574 t, step, wcycle, fplog);
578 rvec *f = as_rvec_array(forceWithVirial->force_.data());
580 /* Add the forces from enforced rotation potentials (if any) */
583 wallcycle_start(wcycle, ewcROTadd);
584 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
585 wallcycle_stop(wcycle, ewcROTadd);
590 /* Note that since init_edsam() is called after the initialization
591 * of forcerec, edsam doesn't request the noVirSum force buffer.
592 * Thus if no other algorithm (e.g. PME) requires it, the forces
593 * here will contribute to the virial.
595 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, bNS);
598 /* Add forces from interactive molecular dynamics (IMD), if any */
599 if (inputrec->bIMD && computeForces)
601 imdSession->applyForces(f);
605 /*! \brief Launch the prepare_step and spread stages of PME GPU.
607 * \param[in] pmedata The PME structure
608 * \param[in] box The box matrix
609 * \param[in] x Coordinate array
610 * \param[in] flags Force flags
611 * \param[in] pmeFlags PME flags
612 * \param[in] wcycle The wallcycle structure
614 static inline void launchPmeGpuSpread(gmx_pme_t *pmedata,
619 gmx_wallcycle_t wcycle)
621 pme_gpu_prepare_computation(pmedata, (flags & GMX_FORCE_DYNAMICBOX) != 0, box, wcycle, pmeFlags);
622 pme_gpu_launch_spread(pmedata, x, wcycle);
625 /*! \brief Launch the FFT and gather stages of PME GPU
627 * This function only implements setting the output forces (no accumulation).
629 * \param[in] pmedata The PME structure
630 * \param[in] wcycle The wallcycle structure
632 static void launchPmeGpuFftAndGather(gmx_pme_t *pmedata,
633 gmx_wallcycle_t wcycle)
635 pme_gpu_launch_complex_transforms(pmedata, wcycle);
636 pme_gpu_launch_gather(pmedata, wcycle, PmeForceOutputHandling::Set);
640 * Polling wait for either of the PME or nonbonded GPU tasks.
642 * Instead of a static order in waiting for GPU tasks, this function
643 * polls checking which of the two tasks completes first, and does the
644 * associated force buffer reduction overlapped with the other task.
645 * By doing that, unlike static scheduling order, it can always overlap
646 * one of the reductions, regardless of the GPU task completion order.
648 * \param[in] nbv Nonbonded verlet structure
649 * \param[in,out] pmedata PME module data
650 * \param[in,out] force Force array to reduce task outputs into.
651 * \param[in,out] forceWithVirial Force and virial buffers
652 * \param[in,out] fshift Shift force output vector results are reduced into
653 * \param[in,out] enerd Energy data structure results are reduced into
654 * \param[in] flags Force flags
655 * \param[in] pmeFlags PME flags
656 * \param[in] haveOtherWork Tells whether there is other work than non-bonded in the stream(s)
657 * \param[in] wcycle The wallcycle structure
659 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t *nbv,
661 gmx::ArrayRefWithPadding<gmx::RVec> *force,
662 gmx::ForceWithVirial *forceWithVirial,
664 gmx_enerdata_t *enerd,
668 gmx_wallcycle_t wcycle)
670 bool isPmeGpuDone = false;
671 bool isNbGpuDone = false;
674 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
676 while (!isPmeGpuDone || !isNbGpuDone)
680 GpuTaskCompletion completionType = (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
681 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, pmeFlags, wcycle, forceWithVirial, enerd, completionType);
686 GpuTaskCompletion completionType = (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
687 wallcycle_start_nocount(wcycle, ewcWAIT_GPU_NB_L);
688 isNbGpuDone = Nbnxm::gpu_try_finish_task(nbv->gpu_nbv,
690 Nbnxm::AtomLocality::Local,
692 enerd->grpp.ener[egLJSR].data(),
693 enerd->grpp.ener[egCOULSR].data(),
694 fshift, completionType);
695 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
696 // To get the call count right, when the task finished we
697 // issue a start/stop.
698 // TODO: move the ewcWAIT_GPU_NB_L cycle counting into nbnxn_gpu_try_finish_task()
699 // and ewcNB_XF_BUF_OPS counting into nbnxn_atomdata_add_nbat_f_to_f().
702 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
703 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
705 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
706 as_rvec_array(force->unpaddedArrayRef().data()), wcycle);
712 /*! \brief Hack structure with force ouput buffers for do_force for the home atoms for this domain */
716 ForceOutputs(rvec *f, gmx::ForceWithVirial const forceWithVirial) :
718 forceWithVirial(forceWithVirial) {}
720 //! Force output buffer used by legacy modules (without SIMD padding)
722 //! Force with direct virial contribution (if there are any; without SIMD padding)
723 gmx::ForceWithVirial forceWithVirial;
726 /*! \brief Set up the different force buffers; also does clearing.
728 * \param[in] fr force record pointer
729 * \param[in] pull_work The pull work object.
730 * \param[in] inputrec input record
731 * \param[in] force force array
732 * \param[in] bDoForces True if force are computed this step
733 * \param[in] doVirial True if virial is computed this step
734 * \param[out] wcycle wallcycle recording structure
736 * \returns Cleared force output structure
739 setupForceOutputs(const t_forcerec *fr,
741 const t_inputrec &inputrec,
742 gmx::ArrayRefWithPadding<gmx::RVec> force,
743 const bool bDoForces,
745 gmx_wallcycle_t wcycle)
747 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
749 /* Temporary solution until all routines take PaddedRVecVector */
750 rvec *const f = as_rvec_array(force.unpaddedArrayRef().data());
753 /* Clear the short- and long-range forces */
754 clear_rvecs_omp(fr->natoms_force_constr, f);
757 /* If we need to compute the virial, we might need a separate
758 * force buffer for algorithms for which the virial is calculated
759 * directly, such as PME. Otherwise, forceWithVirial uses the
760 * the same force (f in legacy calls) buffer as other algorithms.
762 const bool useSeparateForceWithVirialBuffer = (bDoForces && (doVirial && fr->haveDirectVirialContributions));
765 /* forceWithVirial uses the local atom range only */
766 gmx::ForceWithVirial forceWithVirial (useSeparateForceWithVirialBuffer ?
767 *fr->forceBufferForDirectVirialContributions : force.unpaddedArrayRef(),
770 if (useSeparateForceWithVirialBuffer)
772 /* TODO: update comment
773 * We only compute forces on local atoms. Note that vsites can
774 * spread to non-local atoms, but that part of the buffer is
775 * cleared separately in the vsite spreading code.
777 clear_rvecs_omp(forceWithVirial.force_.size(), as_rvec_array(forceWithVirial.force_.data()));
780 if (inputrec.bPull && pull_have_constraint(pull_work))
782 clear_pull_forces(pull_work);
785 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
787 return ForceOutputs(f, forceWithVirial);
791 /*! \brief Set up flags that indicate what type of work is there to compute.
793 * Currently we only update it at search steps,
794 * but some properties may change more frequently (e.g. virial/non-virial step),
795 * so when including those either the frequency of update (per-step) or the scope
796 * of a flag will change (i.e. a set of flags for nstlist steps).
800 setupForceWorkload(gmx::PpForceWorkload *forceWork,
801 const t_inputrec *inputrec,
802 const t_forcerec *fr,
803 const pull_t *pull_work,
810 forceWork->haveSpecialForces = haveSpecialForces(inputrec, fr->forceProviders, pull_work, forceFlags, ed);
811 forceWork->haveCpuBondedWork = haveCpuBondeds(*fr);
812 forceWork->haveGpuBondedWork = ((fr->gpuBonded != nullptr) && fr->gpuBonded->haveInteractions());
813 forceWork->haveRestraintsWork = havePositionRestraints(idef, *fcd);
814 forceWork->haveCpuListedForceWork = haveCpuListedForces(*fr, idef, *fcd);
817 void do_force(FILE *fplog,
819 const gmx_multisim_t *ms,
820 const t_inputrec *inputrec,
822 gmx_enfrot *enforcedRotation,
823 gmx::ImdSession *imdSession,
827 gmx_wallcycle_t wcycle,
828 const gmx_localtop_t *top,
830 gmx::ArrayRefWithPadding<gmx::RVec> x, //NOLINT(performance-unnecessary-value-param)
832 gmx::ArrayRefWithPadding<gmx::RVec> force, //NOLINT(performance-unnecessary-value-param)
834 const t_mdatoms *mdatoms,
835 gmx_enerdata_t *enerd,
837 gmx::ArrayRef<real> lambda,
840 gmx::PpForceWorkload *ppForceWorkload,
841 const gmx_vsite_t *vsite,
846 const DDBalanceRegionHandler &ddBalanceRegionHandler)
850 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
851 gmx_bool bDoForces, bUseGPU, bUseOrEmulGPU;
852 nonbonded_verlet_t *nbv = fr->nbv.get();
853 interaction_const_t *ic = fr->ic;
855 /* modify force flag if not doing nonbonded */
858 flags &= ~GMX_FORCE_NONBONDED;
860 bStateChanged = ((flags & GMX_FORCE_STATECHANGED) != 0);
861 bNS = ((flags & GMX_FORCE_NS) != 0);
862 bFillGrid = (bNS && bStateChanged);
863 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
864 bDoForces = ((flags & GMX_FORCE_FORCES) != 0);
865 bUseGPU = fr->nbv->useGpu();
866 bUseOrEmulGPU = bUseGPU || fr->nbv->emulateGpu();
868 const auto pmeRunMode = fr->pmedata ? pme_run_mode(fr->pmedata) : PmeRunMode::CPU;
869 // TODO slim this conditional down - inputrec and duty checks should mean the same in proper code!
870 const bool useGpuPme = EEL_PME(fr->ic->eeltype) && thisRankHasDuty(cr, DUTY_PME) &&
871 ((pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed));
872 const int pmeFlags = GMX_PME_SPREAD | GMX_PME_SOLVE |
873 ((flags & GMX_FORCE_VIRIAL) ? GMX_PME_CALC_ENER_VIR : 0) |
874 ((flags & GMX_FORCE_ENERGY) ? GMX_PME_CALC_ENER_VIR : 0) |
875 ((flags & GMX_FORCE_FORCES) ? GMX_PME_CALC_F : 0);
877 const bool useGpuXBufOps = (c_enableGpuBufOps && bUseGPU && (GMX_GPU == GMX_GPU_CUDA));
879 /* At a search step we need to start the first balancing region
880 * somewhere early inside the step after communication during domain
881 * decomposition (and not during the previous step as usual).
885 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
889 const int homenr = mdatoms->homenr;
891 clear_mat(vir_force);
895 update_forcerec(fr, box);
897 if (inputrecNeedMutot(inputrec))
899 /* Calculate total (local) dipole moment in a temporary common array.
900 * This makes it possible to sum them over nodes faster.
902 calc_mu(start, homenr,
903 x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
908 if (fr->ePBC != epbcNONE)
910 /* Compute shift vectors every step,
911 * because of pressure coupling or box deformation!
913 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
915 calc_shifts(box, fr->shift_vec);
920 put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr), gmx_omp_nthreads_get(emntDefault));
921 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
923 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
925 unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
929 nbnxn_atomdata_copy_shiftvec((flags & GMX_FORCE_DYNAMICBOX) != 0,
930 fr->shift_vec, nbv->nbat.get());
933 if (!thisRankHasDuty(cr, DUTY_PME))
935 /* Send particle coordinates to the pme nodes.
936 * Since this is only implemented for domain decomposition
937 * and domain decomposition does not use the graph,
938 * we do not need to worry about shifting.
940 gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
941 lambda[efptCOUL], lambda[efptVDW],
942 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)) != 0,
949 launchPmeGpuSpread(fr->pmedata, box, as_rvec_array(x.unpaddedArrayRef().data()), flags, pmeFlags, wcycle);
952 /* do gridding for pair search */
955 if (graph && bStateChanged)
957 /* Calculate intramolecular shift vectors to make molecules whole */
958 mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
962 // - vzero is constant, do we need to pass it?
963 // - box_diag should be passed directly to nbnxn_put_on_grid
969 box_diag[XX] = box[XX][XX];
970 box_diag[YY] = box[YY][YY];
971 box_diag[ZZ] = box[ZZ][ZZ];
973 wallcycle_start(wcycle, ewcNS);
974 if (!DOMAINDECOMP(cr))
976 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
977 nbnxn_put_on_grid(nbv, box,
979 nullptr, 0, mdatoms->homenr, -1,
980 fr->cginfo, x.unpaddedArrayRef(),
982 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
986 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
987 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd),
988 fr->cginfo, x.unpaddedArrayRef());
989 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
992 nbv->setAtomProperties(*mdatoms, fr->cginfo);
994 wallcycle_stop(wcycle, ewcNS);
996 /* initialize the GPU nbnxm atom data and bonded data structures */
999 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1001 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1002 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1003 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1007 /* Now we put all atoms on the grid, we can assign bonded
1008 * interactions to the GPU, where the grid order is
1009 * needed. Also the xq, f and fshift device buffers have
1010 * been reallocated if needed, so the bonded code can
1011 * learn about them. */
1012 // TODO the xq, f, and fshift buffers are now shared
1013 // resources, so they should be maintained by a
1014 // higher-level object than the nb module.
1015 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1017 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1018 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1019 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1021 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1024 // Need to run after the GPU-offload bonded interaction lists
1025 // are set up to be able to determine whether there is bonded work.
1026 setupForceWorkload(ppForceWorkload,
1036 /* do local pair search */
1039 // TODO: fuse this branch with the above bNS block
1040 wallcycle_start_nocount(wcycle, ewcNS);
1041 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1042 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1043 nbv->constructPairlist(Nbnxm::InteractionLocality::Local,
1044 &top->excls, step, nrnb);
1045 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1046 wallcycle_stop(wcycle, ewcNS);
1050 nbv->atomdata_init_copy_x_to_nbat_x_gpu( Nbnxm::AtomLocality::Local);
1056 nbv->setCoordinates(Nbnxm::AtomLocality::Local, false,
1057 x.unpaddedArrayRef(), useGpuXBufOps, pme_gpu_get_device_x(fr->pmedata), wcycle);
1062 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1064 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1066 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1067 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1068 if (bNS || !useGpuXBufOps)
1070 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1071 Nbnxm::AtomLocality::Local,
1072 ppForceWorkload->haveGpuBondedWork);
1074 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1075 // with X buffer ops offloaded to the GPU on all but the search steps
1077 // bonded work not split into separate local and non-local, so with DD
1078 // we can only launch the kernel after non-local coordinates have been received.
1079 if (ppForceWorkload->haveGpuBondedWork && !havePPDomainDecomposition(cr))
1081 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1082 fr->gpuBonded->launchKernel(fr, flags, box);
1083 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1086 /* launch local nonbonded work on GPU */
1087 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1088 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFNo,
1089 step, nrnb, wcycle);
1090 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1091 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1096 // In PME GPU and mixed mode we launch FFT / gather after the
1097 // X copy/transform to allow overlap as well as after the GPU NB
1098 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1099 // the nonbonded kernel.
1100 launchPmeGpuFftAndGather(fr->pmedata, wcycle);
1103 /* Communicate coordinates and sum dipole if necessary +
1104 do non-local pair search */
1105 if (havePPDomainDecomposition(cr))
1109 // TODO: fuse this branch with the above large bNS block
1110 wallcycle_start_nocount(wcycle, ewcNS);
1111 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1112 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1113 nbv->constructPairlist(Nbnxm::InteractionLocality::NonLocal,
1114 &top->excls, step, nrnb);
1115 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1116 wallcycle_stop(wcycle, ewcNS);
1121 nbv->atomdata_init_copy_x_to_nbat_x_gpu( Nbnxm::AtomLocality::NonLocal);
1126 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1128 nbv->setCoordinates(Nbnxm::AtomLocality::NonLocal, false,
1129 x.unpaddedArrayRef(), useGpuXBufOps, pme_gpu_get_device_x(fr->pmedata), wcycle);
1135 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1137 if (bNS || !useGpuXBufOps)
1139 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1140 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1141 Nbnxm::AtomLocality::NonLocal,
1142 ppForceWorkload->haveGpuBondedWork);
1143 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1146 if (ppForceWorkload->haveGpuBondedWork)
1148 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1149 fr->gpuBonded->launchKernel(fr, flags, box);
1150 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1153 /* launch non-local nonbonded tasks on GPU */
1154 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1155 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1156 step, nrnb, wcycle);
1157 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1159 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1165 /* launch D2H copy-back F */
1166 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1167 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1168 if (havePPDomainDecomposition(cr))
1170 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1171 flags, Nbnxm::AtomLocality::NonLocal, ppForceWorkload->haveGpuBondedWork);
1173 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1174 flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork);
1175 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1177 if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1179 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1180 fr->gpuBonded->launchEnergyTransfer();
1181 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1183 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1186 if (bStateChanged && inputrecNeedMutot(inputrec))
1190 gmx_sumd(2*DIM, mu, cr);
1192 ddBalanceRegionHandler.reopenRegionCpu();
1195 for (i = 0; i < 2; i++)
1197 for (j = 0; j < DIM; j++)
1199 fr->mu_tot[i][j] = mu[i*DIM + j];
1203 if (fr->efep == efepNO)
1205 copy_rvec(fr->mu_tot[0], mu_tot);
1209 for (j = 0; j < DIM; j++)
1212 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1213 lambda[efptCOUL]*fr->mu_tot[1][j];
1217 /* Reset energies */
1218 reset_enerdata(enerd);
1219 clear_rvecs(SHIFTS, fr->fshift);
1221 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1223 wallcycle_start(wcycle, ewcPPDURINGPME);
1224 dd_force_flop_start(cr->dd, nrnb);
1229 wallcycle_start(wcycle, ewcROT);
1230 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, bNS);
1231 wallcycle_stop(wcycle, ewcROT);
1234 /* Start the force cycle counter.
1235 * Note that a different counter is used for dynamic load balancing.
1237 wallcycle_start(wcycle, ewcFORCE);
1239 // set up and clear force outputs
1240 struct ForceOutputs forceOut = setupForceOutputs(fr, pull_work, *inputrec, force, bDoForces,
1241 ((flags & GMX_FORCE_VIRIAL) != 0), wcycle);
1243 /* We calculate the non-bonded forces, when done on the CPU, here.
1244 * We do this before calling do_force_lowlevel, because in that
1245 * function, the listed forces are calculated before PME, which
1246 * does communication. With this order, non-bonded and listed
1247 * force calculation imbalance can be balanced out by the domain
1248 * decomposition load balancing.
1253 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFYes,
1254 step, nrnb, wcycle);
1257 if (fr->efep != efepNO)
1259 /* Calculate the local and non-local free energy interactions here.
1260 * Happens here on the CPU both with and without GPU.
1262 wallcycle_sub_start(wcycle, ewcsNONBONDED);
1263 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::Local,
1264 fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1265 inputrec->fepvals, lambda.data(),
1266 enerd, flags, nrnb);
1268 if (havePPDomainDecomposition(cr))
1270 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::NonLocal,
1271 fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1272 inputrec->fepvals, lambda.data(),
1273 enerd, flags, nrnb);
1275 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
1280 if (havePPDomainDecomposition(cr))
1282 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1283 step, nrnb, wcycle);
1286 /* Add all the non-bonded force to the normal force array.
1287 * This can be split into a local and a non-local part when overlapping
1288 * communication with calculation with domain decomposition.
1290 wallcycle_stop(wcycle, ewcFORCE);
1292 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::All, forceOut.f, wcycle);
1294 wallcycle_start_nocount(wcycle, ewcFORCE);
1296 /* If there are multiple fshift output buffers we need to reduce them */
1297 if (flags & GMX_FORCE_VIRIAL)
1299 /* This is not in a subcounter because it takes a
1300 negligible and constant-sized amount of time */
1301 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->nbat.get(),
1306 /* update QMMMrec, if necessary */
1309 update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1312 /* Compute the bonded and non-bonded energies and optionally forces */
1313 do_force_lowlevel(fr, inputrec, &(top->idef),
1314 cr, ms, nrnb, wcycle, mdatoms,
1315 x, hist, forceOut.f, &forceOut.forceWithVirial, enerd, fcd,
1316 box, lambda.data(), graph, fr->mu_tot,
1318 ddBalanceRegionHandler);
1320 wallcycle_stop(wcycle, ewcFORCE);
1322 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1323 imdSession, pull_work, step, t, wcycle,
1324 fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda.data(),
1325 flags, &forceOut.forceWithVirial, enerd,
1328 // Will store the amount of cycles spent waiting for the GPU that
1329 // will be later used in the DLB accounting.
1330 float cycles_wait_gpu = 0;
1333 /* wait for non-local forces (or calculate in emulation mode) */
1334 if (havePPDomainDecomposition(cr))
1338 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1339 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1340 flags, Nbnxm::AtomLocality::NonLocal,
1341 ppForceWorkload->haveGpuBondedWork,
1342 enerd->grpp.ener[egLJSR].data(),
1343 enerd->grpp.ener[egCOULSR].data(),
1345 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1349 wallcycle_start_nocount(wcycle, ewcFORCE);
1350 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFYes,
1351 step, nrnb, wcycle);
1352 wallcycle_stop(wcycle, ewcFORCE);
1355 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::NonLocal,
1356 forceOut.f, wcycle);
1358 if (fr->nbv->emulateGpu() && (flags & GMX_FORCE_VIRIAL))
1360 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->nbat.get(),
1366 if (havePPDomainDecomposition(cr))
1368 /* We are done with the CPU compute.
1369 * We will now communicate the non-local forces.
1370 * If we use a GPU this will overlap with GPU work, so in that case
1371 * we do not close the DD force balancing region here.
1373 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1377 dd_move_f(cr->dd, force.unpaddedArrayRef(), fr->fshift, wcycle);
1381 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1382 // an alternating wait/reduction scheme.
1383 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPme && bUseGPU && !DOMAINDECOMP(cr));
1384 if (alternateGpuWait)
1386 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &force, &forceOut.forceWithVirial, fr->fshift, enerd,
1387 flags, pmeFlags, ppForceWorkload->haveGpuBondedWork, wcycle);
1390 if (!alternateGpuWait && useGpuPme)
1392 pme_gpu_wait_and_reduce(fr->pmedata, pmeFlags, wcycle, &forceOut.forceWithVirial, enerd);
1395 /* Wait for local GPU NB outputs on the non-alternating wait path */
1396 if (!alternateGpuWait && bUseGPU)
1398 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1399 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1400 * but even with a step of 0.1 ms the difference is less than 1%
1403 const float gpuWaitApiOverheadMargin = 2e6f; /* cycles */
1405 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1406 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1407 flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork,
1408 enerd->grpp.ener[egLJSR].data(),
1409 enerd->grpp.ener[egCOULSR].data(),
1411 float cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1413 if (ddBalanceRegionHandler.useBalancingRegion())
1415 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1416 if (bDoForces && cycles_tmp <= gpuWaitApiOverheadMargin)
1418 /* We measured few cycles, it could be that the kernel
1419 * and transfer finished earlier and there was no actual
1420 * wait time, only API call overhead.
1421 * Then the actual time could be anywhere between 0 and
1422 * cycles_wait_est. We will use half of cycles_wait_est.
1424 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1426 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1430 if (fr->nbv->emulateGpu())
1432 // NOTE: emulation kernel is not included in the balancing region,
1433 // but emulation mode does not target performance anyway
1434 wallcycle_start_nocount(wcycle, ewcFORCE);
1435 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local,
1436 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1437 step, nrnb, wcycle);
1438 wallcycle_stop(wcycle, ewcFORCE);
1443 pme_gpu_reinit_computation(fr->pmedata, wcycle);
1448 /* now clear the GPU outputs while we finish the step on the CPU */
1449 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1450 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1451 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, flags);
1453 if (nbv->isDynamicPruningStepGpu(step))
1455 nbv->dispatchPruneKernelGpu(step);
1457 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1458 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1461 if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1463 wallcycle_start(wcycle, ewcWAIT_GPU_BONDED);
1464 // in principle this should be included in the DD balancing region,
1465 // but generally it is infrequent so we'll omit it for the sake of
1467 fr->gpuBonded->accumulateEnergyTerms(enerd);
1468 wallcycle_stop(wcycle, ewcWAIT_GPU_BONDED);
1470 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1471 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1472 fr->gpuBonded->clearEnergies();
1473 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1474 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1477 /* Do the nonbonded GPU (or emulation) force buffer reduction
1478 * on the non-alternating path. */
1479 if (bUseOrEmulGPU && !alternateGpuWait)
1481 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
1482 forceOut.f, wcycle);
1484 if (DOMAINDECOMP(cr))
1486 dd_force_flop_stop(cr->dd, nrnb);
1491 /* If we have NoVirSum forces, but we do not calculate the virial,
1492 * we sum fr->f_novirsum=forceOut.f later.
1494 if (vsite && !(fr->haveDirectVirialContributions && !(flags & GMX_FORCE_VIRIAL)))
1496 spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, fr->fshift, FALSE, nullptr, nrnb,
1497 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1500 if (flags & GMX_FORCE_VIRIAL)
1502 /* Calculation of the virial must be done after vsites! */
1503 calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f,
1504 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1508 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1510 /* In case of node-splitting, the PP nodes receive the long-range
1511 * forces, virial and energy from the PME nodes here.
1513 pme_receive_force_ener(cr, &forceOut.forceWithVirial, enerd, wcycle);
1518 post_process_forces(cr, step, nrnb, wcycle,
1519 top, box, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, &forceOut.forceWithVirial,
1520 vir_force, mdatoms, graph, fr, vsite,
1524 if (flags & GMX_FORCE_ENERGY)
1526 /* Sum the potential energy terms from group contributions */
1527 sum_epot(&(enerd->grpp), enerd->term);
1529 if (!EI_TPI(inputrec->eI))
1531 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1535 /* In case we don't have constraints and are using GPUs, the next balancing
1536 * region starts here.
1537 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1538 * virial calculation and COM pulling, is not thus not included in
1539 * the balance timing, which is ok as most tasks do communication.
1541 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);