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 using gmx::ForceOutputs;
119 // TODO: this environment variable allows us to verify before release
120 // that on less common architectures the total cost of polling is not larger than
121 // a blocking wait (so polling does not introduce overhead when the static
122 // PME-first ordering would suffice).
123 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
125 // environment variable to enable GPU buffer ops, to allow incremental and optional
126 // introduction of this functionality.
127 // TODO eventially tie this in with other existing GPU flags.
128 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,
262 ForceOutputs *forceOutputs,
264 const t_mdatoms *mdatoms,
265 const t_graph *graph,
266 const t_forcerec *fr,
267 const gmx_vsite_t *vsite,
270 rvec *f = forceOutputs->f();
272 if (fr->haveDirectVirialContributions)
274 auto &forceWithVirial = forceOutputs->forceWithVirial();
275 rvec *fDirectVir = as_rvec_array(forceWithVirial.force_.data());
279 /* Spread the mesh force on virtual sites to the other particles...
280 * This is parallellized. MPI communication is performed
281 * if the constructing atoms aren't local.
283 matrix virial = { { 0 } };
284 spread_vsite_f(vsite, x, fDirectVir, nullptr,
285 (flags & GMX_FORCE_VIRIAL) != 0, virial,
287 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
288 forceWithVirial.addVirialContribution(virial);
291 if (flags & GMX_FORCE_VIRIAL)
293 /* Now add the forces, this is local */
294 sum_forces(f, forceWithVirial.force_);
296 /* Add the direct virial contributions */
297 GMX_ASSERT(forceWithVirial.computeVirial_, "forceWithVirial should request virial computation when we request the virial");
298 m_add(vir_force, forceWithVirial.getVirial(), vir_force);
302 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
307 if (fr->print_force >= 0)
309 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
313 static void do_nb_verlet(t_forcerec *fr,
314 const interaction_const_t *ic,
315 gmx_enerdata_t *enerd,
317 const Nbnxm::InteractionLocality ilocality,
321 gmx_wallcycle_t wcycle)
323 if (!(flags & GMX_FORCE_NONBONDED))
325 /* skip non-bonded calculation */
329 nonbonded_verlet_t *nbv = fr->nbv.get();
331 /* GPU kernel launch overhead is already timed separately */
332 if (fr->cutoff_scheme != ecutsVERLET)
334 gmx_incons("Invalid cut-off scheme passed!");
339 /* When dynamic pair-list pruning is requested, we need to prune
340 * at nstlistPrune steps.
342 if (nbv->isDynamicPruningStepCpu(step))
344 /* Prune the pair-list beyond fr->ic->rlistPrune using
345 * the current coordinates of the atoms.
347 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
348 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
349 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
353 nbv->dispatchNonbondedKernel(ilocality, *ic, flags, clearF, *fr, enerd, nrnb);
356 static inline void clear_rvecs_omp(int n, rvec v[])
358 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, n);
360 /* Note that we would like to avoid this conditional by putting it
361 * into the omp pragma instead, but then we still take the full
362 * omp parallel for overhead (at least with gcc5).
366 for (int i = 0; i < n; i++)
373 #pragma omp parallel for num_threads(nth) schedule(static)
374 for (int i = 0; i < n; i++)
381 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
383 * \param groupOptions Group options, containing T-coupling options
385 static real averageKineticEnergyEstimate(const t_grpopts &groupOptions)
387 real nrdfCoupled = 0;
388 real nrdfUncoupled = 0;
389 real kineticEnergy = 0;
390 for (int g = 0; g < groupOptions.ngtc; g++)
392 if (groupOptions.tau_t[g] >= 0)
394 nrdfCoupled += groupOptions.nrdf[g];
395 kineticEnergy += groupOptions.nrdf[g]*0.5*groupOptions.ref_t[g]*BOLTZ;
399 nrdfUncoupled += groupOptions.nrdf[g];
403 /* This conditional with > also catches nrdf=0 */
404 if (nrdfCoupled > nrdfUncoupled)
406 return kineticEnergy*(nrdfCoupled + nrdfUncoupled)/nrdfCoupled;
414 /*! \brief This routine checks that the potential energy is finite.
416 * Always checks that the potential energy is finite. If step equals
417 * inputrec.init_step also checks that the magnitude of the potential energy
418 * is reasonable. Terminates with a fatal error when a check fails.
419 * Note that passing this check does not guarantee finite forces,
420 * since those use slightly different arithmetics. But in most cases
421 * there is just a narrow coordinate range where forces are not finite
422 * and energies are finite.
424 * \param[in] step The step number, used for checking and printing
425 * \param[in] enerd The energy data; the non-bonded group energies need to be added to enerd.term[F_EPOT] before calling this routine
426 * \param[in] inputrec The input record
428 static void checkPotentialEnergyValidity(int64_t step,
429 const gmx_enerdata_t &enerd,
430 const t_inputrec &inputrec)
432 /* Threshold valid for comparing absolute potential energy against
433 * the kinetic energy. Normally one should not consider absolute
434 * potential energy values, but with a factor of one million
435 * we should never get false positives.
437 constexpr real c_thresholdFactor = 1e6;
439 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
440 real averageKineticEnergy = 0;
441 /* We only check for large potential energy at the initial step,
442 * because that is by far the most likely step for this too occur
443 * and because computing the average kinetic energy is not free.
444 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
445 * before they become NaN.
447 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
449 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
452 if (energyIsNotFinite || (averageKineticEnergy > 0 &&
453 enerd.term[F_EPOT] > c_thresholdFactor*averageKineticEnergy))
455 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.",
458 energyIsNotFinite ? "not finite" : "extremely high",
460 enerd.term[F_COUL_SR],
461 energyIsNotFinite ? "non-finite" : "very high",
462 energyIsNotFinite ? " or Nan" : "");
466 /*! \brief Return true if there are special forces computed this step.
468 * The conditionals exactly correspond to those in computeSpecialForces().
471 haveSpecialForces(const t_inputrec *inputrec,
472 ForceProviders *forceProviders,
473 const pull_t *pull_work,
477 const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
480 ((computeForces && forceProviders->hasForceProvider()) || // forceProviders
481 (inputrec->bPull && pull_have_potential(pull_work)) || // pull
482 inputrec->bRot || // enforced rotation
483 (ed != nullptr) || // flooding
484 (inputrec->bIMD && computeForces)); // IMD
487 /*! \brief Compute forces and/or energies for special algorithms
489 * The intention is to collect all calls to algorithms that compute
490 * forces on local atoms only and that do not contribute to the local
491 * virial sum (but add their virial contribution separately).
492 * Eventually these should likely all become ForceProviders.
493 * Within this function the intention is to have algorithms that do
494 * global communication at the end, so global barriers within the MD loop
495 * are as close together as possible.
497 * \param[in] fplog The log file
498 * \param[in] cr The communication record
499 * \param[in] inputrec The input record
500 * \param[in] awh The Awh module (nullptr if none in use).
501 * \param[in] enforcedRotation Enforced rotation module.
502 * \param[in] imdSession The IMD session
503 * \param[in] pull_work The pull work structure.
504 * \param[in] step The current MD step
505 * \param[in] t The current time
506 * \param[in,out] wcycle Wallcycle accounting struct
507 * \param[in,out] forceProviders Pointer to a list of force providers
508 * \param[in] box The unit cell
509 * \param[in] x The coordinates
510 * \param[in] mdatoms Per atom properties
511 * \param[in] lambda Array of free-energy lambda values
512 * \param[in] forceFlags Flags that tell whether we should compute forces/energies/virial
513 * \param[in,out] forceWithVirial Force and virial buffers
514 * \param[in,out] enerd Energy buffer
515 * \param[in,out] ed Essential dynamics pointer
516 * \param[in] bNS Tells if we did neighbor searching this step, used for ED sampling
518 * \todo Remove bNS, which is used incorrectly.
519 * \todo Convert all other algorithms called here to ForceProviders.
522 computeSpecialForces(FILE *fplog,
524 const t_inputrec *inputrec,
526 gmx_enfrot *enforcedRotation,
527 gmx::ImdSession *imdSession,
531 gmx_wallcycle_t wcycle,
532 ForceProviders *forceProviders,
534 gmx::ArrayRef<const gmx::RVec> x,
535 const t_mdatoms *mdatoms,
538 gmx::ForceWithVirial *forceWithVirial,
539 gmx_enerdata_t *enerd,
543 const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
545 /* NOTE: Currently all ForceProviders only provide forces.
546 * When they also provide energies, remove this conditional.
550 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
551 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
553 /* Collect forces from modules */
554 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
557 if (inputrec->bPull && pull_have_potential(pull_work))
559 pull_potential_wrapper(cr, inputrec, box, x,
561 mdatoms, enerd, pull_work, lambda, t,
566 enerd->term[F_COM_PULL] +=
567 awh->applyBiasForcesAndUpdateBias(inputrec->ePBC, *mdatoms, box,
569 t, step, wcycle, fplog);
573 rvec *f = as_rvec_array(forceWithVirial->force_.data());
575 /* Add the forces from enforced rotation potentials (if any) */
578 wallcycle_start(wcycle, ewcROTadd);
579 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
580 wallcycle_stop(wcycle, ewcROTadd);
585 /* Note that since init_edsam() is called after the initialization
586 * of forcerec, edsam doesn't request the noVirSum force buffer.
587 * Thus if no other algorithm (e.g. PME) requires it, the forces
588 * here will contribute to the virial.
590 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, bNS);
593 /* Add forces from interactive molecular dynamics (IMD), if any */
594 if (inputrec->bIMD && computeForces)
596 imdSession->applyForces(f);
600 /*! \brief Launch the prepare_step and spread stages of PME GPU.
602 * \param[in] pmedata The PME structure
603 * \param[in] box The box matrix
604 * \param[in] x Coordinate array
605 * \param[in] flags Force flags
606 * \param[in] pmeFlags PME flags
607 * \param[in] wcycle The wallcycle structure
609 static inline void launchPmeGpuSpread(gmx_pme_t *pmedata,
614 gmx_wallcycle_t wcycle)
616 pme_gpu_prepare_computation(pmedata, (flags & GMX_FORCE_DYNAMICBOX) != 0, box, wcycle, pmeFlags);
617 pme_gpu_launch_spread(pmedata, x, wcycle);
620 /*! \brief Launch the FFT and gather stages of PME GPU
622 * This function only implements setting the output forces (no accumulation).
624 * \param[in] pmedata The PME structure
625 * \param[in] wcycle The wallcycle structure
626 * \param[in] useGpuFPmeReduction Whether forces will be reduced on GPU
628 static void launchPmeGpuFftAndGather(gmx_pme_t *pmedata,
629 gmx_wallcycle_t wcycle,
630 bool useGpuFPmeReduction)
632 pme_gpu_launch_complex_transforms(pmedata, wcycle);
633 pme_gpu_launch_gather(pmedata, wcycle, PmeForceOutputHandling::Set, useGpuFPmeReduction);
637 * Polling wait for either of the PME or nonbonded GPU tasks.
639 * Instead of a static order in waiting for GPU tasks, this function
640 * polls checking which of the two tasks completes first, and does the
641 * associated force buffer reduction overlapped with the other task.
642 * By doing that, unlike static scheduling order, it can always overlap
643 * one of the reductions, regardless of the GPU task completion order.
645 * \param[in] nbv Nonbonded verlet structure
646 * \param[in,out] pmedata PME module data
647 * \param[in,out] force Force array to reduce task outputs into.
648 * \param[in,out] forceWithVirial Force and virial buffers
649 * \param[in,out] fshift Shift force output vector results are reduced into
650 * \param[in,out] enerd Energy data structure results are reduced into
651 * \param[in] flags Force flags
652 * \param[in] pmeFlags PME flags
653 * \param[in] wcycle The wallcycle structure
655 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t *nbv,
657 gmx::ArrayRefWithPadding<gmx::RVec> *force,
658 gmx::ForceWithVirial *forceWithVirial,
660 gmx_enerdata_t *enerd,
663 gmx_wallcycle_t wcycle)
665 bool isPmeGpuDone = false;
666 bool isNbGpuDone = false;
669 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
671 while (!isPmeGpuDone || !isNbGpuDone)
675 GpuTaskCompletion completionType = (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
676 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, pmeFlags, wcycle, forceWithVirial, enerd, completionType);
681 GpuTaskCompletion completionType = (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
682 wallcycle_start_nocount(wcycle, ewcWAIT_GPU_NB_L);
683 isNbGpuDone = Nbnxm::gpu_try_finish_task(nbv->gpu_nbv,
685 Nbnxm::AtomLocality::Local,
686 enerd->grpp.ener[egLJSR].data(),
687 enerd->grpp.ener[egCOULSR].data(),
688 fshift, completionType);
689 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
690 // To get the call count right, when the task finished we
691 // issue a start/stop.
692 // TODO: move the ewcWAIT_GPU_NB_L cycle counting into nbnxn_gpu_try_finish_task()
693 // and ewcNB_XF_BUF_OPS counting into nbnxn_atomdata_add_nbat_f_to_f().
696 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
697 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
699 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
700 as_rvec_array(force->unpaddedArrayRef().data()));
706 /*! \brief Set up the different force buffers; also does clearing.
708 * \param[in] fr force record pointer
709 * \param[in] pull_work The pull work object.
710 * \param[in] inputrec input record
711 * \param[in] force force array
712 * \param[in] bDoForces True if force are computed this step
713 * \param[in] doVirial True if virial is computed this step
714 * \param[out] wcycle wallcycle recording structure
716 * \returns Cleared force output structure
719 setupForceOutputs(t_forcerec *fr,
721 const t_inputrec &inputrec,
722 gmx::ArrayRefWithPadding<gmx::RVec> force,
723 const bool bDoForces,
725 gmx_wallcycle_t wcycle)
727 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
729 // TODO: Add the shift force buffer and use it
730 gmx::ForceWithShiftForces forceWithShiftForces(force, doVirial, gmx::ArrayRef<gmx::RVec>());
734 /* Clear the short- and long-range forces */
735 clear_rvecs_omp(fr->natoms_force_constr, forceWithShiftForces.f());
738 /* If we need to compute the virial, we might need a separate
739 * force buffer for algorithms for which the virial is calculated
740 * directly, such as PME. Otherwise, forceWithVirial uses the
741 * the same force (f in legacy calls) buffer as other algorithms.
743 const bool useSeparateForceWithVirialBuffer = (bDoForces && (doVirial && fr->haveDirectVirialContributions));
745 /* forceWithVirial uses the local atom range only */
746 gmx::ForceWithVirial forceWithVirial(useSeparateForceWithVirialBuffer ?
747 fr->forceBufferForDirectVirialContributions : force.unpaddedArrayRef(),
750 if (useSeparateForceWithVirialBuffer)
752 /* TODO: update comment
753 * We only compute forces on local atoms. Note that vsites can
754 * spread to non-local atoms, but that part of the buffer is
755 * cleared separately in the vsite spreading code.
757 clear_rvecs_omp(forceWithVirial.force_.size(), as_rvec_array(forceWithVirial.force_.data()));
760 if (inputrec.bPull && pull_have_constraint(pull_work))
762 clear_pull_forces(pull_work);
765 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
767 return ForceOutputs(forceWithShiftForces, forceWithVirial);
771 /*! \brief Set up flags that indicate what type of work is there to compute.
773 * Currently we only update it at search steps,
774 * but some properties may change more frequently (e.g. virial/non-virial step),
775 * so when including those either the frequency of update (per-step) or the scope
776 * of a flag will change (i.e. a set of flags for nstlist steps).
780 setupForceWorkload(gmx::PpForceWorkload *forceWork,
781 const t_inputrec *inputrec,
782 const t_forcerec *fr,
783 const pull_t *pull_work,
790 forceWork->haveSpecialForces = haveSpecialForces(inputrec, fr->forceProviders, pull_work, forceFlags, ed);
791 forceWork->haveCpuBondedWork = haveCpuBondeds(*fr);
792 forceWork->haveGpuBondedWork = ((fr->gpuBonded != nullptr) && fr->gpuBonded->haveInteractions());
793 forceWork->haveRestraintsWork = havePositionRestraints(idef, *fcd);
794 forceWork->haveCpuListedForceWork = haveCpuListedForces(*fr, idef, *fcd);
798 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
800 * TODO: eliminate the \p useGpuNonbonded and \p useGpuNonbonded when these are
801 * incorporated in PpForceWorkload.
804 launchGpuEndOfStepTasks(nonbonded_verlet_t *nbv,
805 gmx::GpuBonded *gpuBonded,
807 gmx_enerdata_t *enerd,
808 const gmx::PpForceWorkload &forceWorkload,
809 bool useGpuNonbonded,
813 gmx_wallcycle_t wcycle)
817 /* Launch pruning before buffer clearing because the API overhead of the
818 * clear kernel launches can leave the GPU idle while it could be running
821 if (nbv->isDynamicPruningStepGpu(step))
823 nbv->dispatchPruneKernelGpu(step);
826 /* now clear the GPU outputs while we finish the step on the CPU */
827 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
828 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
829 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, flags);
830 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
831 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
836 pme_gpu_reinit_computation(pmedata, wcycle);
839 if (forceWorkload.haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
841 // in principle this should be included in the DD balancing region,
842 // but generally it is infrequent so we'll omit it for the sake of
844 gpuBonded->waitAccumulateEnergyTerms(enerd);
846 gpuBonded->clearEnergies();
851 void do_force(FILE *fplog,
853 const gmx_multisim_t *ms,
854 const t_inputrec *inputrec,
856 gmx_enfrot *enforcedRotation,
857 gmx::ImdSession *imdSession,
861 gmx_wallcycle_t wcycle,
862 const gmx_localtop_t *top,
864 gmx::ArrayRefWithPadding<gmx::RVec> x, //NOLINT(performance-unnecessary-value-param)
866 gmx::ArrayRefWithPadding<gmx::RVec> force, //NOLINT(performance-unnecessary-value-param)
868 const t_mdatoms *mdatoms,
869 gmx_enerdata_t *enerd,
871 gmx::ArrayRef<real> lambda,
874 gmx::PpForceWorkload *ppForceWorkload,
875 const gmx_vsite_t *vsite,
880 const DDBalanceRegionHandler &ddBalanceRegionHandler)
884 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
885 gmx_bool bDoForces, bUseGPU, bUseOrEmulGPU;
886 nonbonded_verlet_t *nbv = fr->nbv.get();
887 interaction_const_t *ic = fr->ic;
889 /* modify force flag if not doing nonbonded */
892 flags &= ~GMX_FORCE_NONBONDED;
894 bStateChanged = ((flags & GMX_FORCE_STATECHANGED) != 0);
895 bNS = ((flags & GMX_FORCE_NS) != 0);
896 bFillGrid = (bNS && bStateChanged);
897 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
898 bDoForces = ((flags & GMX_FORCE_FORCES) != 0);
899 bUseGPU = fr->nbv->useGpu();
900 bUseOrEmulGPU = bUseGPU || fr->nbv->emulateGpu();
902 const auto pmeRunMode = fr->pmedata ? pme_run_mode(fr->pmedata) : PmeRunMode::CPU;
903 // TODO slim this conditional down - inputrec and duty checks should mean the same in proper code!
904 const bool useGpuPme = EEL_PME(fr->ic->eeltype) && thisRankHasDuty(cr, DUTY_PME) &&
905 ((pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed));
906 const int pmeFlags = GMX_PME_SPREAD | GMX_PME_SOLVE |
907 ((flags & GMX_FORCE_VIRIAL) ? GMX_PME_CALC_ENER_VIR : 0) |
908 ((flags & GMX_FORCE_ENERGY) ? GMX_PME_CALC_ENER_VIR : 0) |
909 ((flags & GMX_FORCE_FORCES) ? GMX_PME_CALC_F : 0);
911 // Switches on whether to use GPU for position and force buffer operations
912 // TODO consider all possible combinations of triggers, and how to combine optimally in each case.
913 const BufferOpsUseGpu useGpuXBufOps = (c_enableGpuBufOps && bUseGPU && (GMX_GPU == GMX_GPU_CUDA)) ?
914 BufferOpsUseGpu::True : BufferOpsUseGpu::False;;
915 // GPU Force buffer ops are disabled on virial steps, because the virial calc is not yet ported to GPU
916 const BufferOpsUseGpu useGpuFBufOps = (c_enableGpuBufOps && bUseGPU && (GMX_GPU == GMX_GPU_CUDA))
917 && !(flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)) ?
918 BufferOpsUseGpu::True : BufferOpsUseGpu::False;
919 // TODO: move / add this flag to the internal PME GPU data structures
920 const bool useGpuFPmeReduction = (useGpuFBufOps == BufferOpsUseGpu::True) &&
921 thisRankHasDuty(cr, DUTY_PME) && useGpuPme; // only supported if this rank is perfoming PME on the GPU
923 /* At a search step we need to start the first balancing region
924 * somewhere early inside the step after communication during domain
925 * decomposition (and not during the previous step as usual).
929 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
933 const int homenr = mdatoms->homenr;
935 clear_mat(vir_force);
939 if (inputrecNeedMutot(inputrec))
941 /* Calculate total (local) dipole moment in a temporary common array.
942 * This makes it possible to sum them over nodes faster.
944 calc_mu(start, homenr,
945 x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
950 if (fr->ePBC != epbcNONE)
952 /* Compute shift vectors every step,
953 * because of pressure coupling or box deformation!
955 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
957 calc_shifts(box, fr->shift_vec);
962 put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr), gmx_omp_nthreads_get(emntDefault));
963 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
965 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
967 unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
971 nbnxn_atomdata_copy_shiftvec((flags & GMX_FORCE_DYNAMICBOX) != 0,
972 fr->shift_vec, nbv->nbat.get());
975 if (!thisRankHasDuty(cr, DUTY_PME))
977 /* Send particle coordinates to the pme nodes.
978 * Since this is only implemented for domain decomposition
979 * and domain decomposition does not use the graph,
980 * we do not need to worry about shifting.
982 gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
983 lambda[efptCOUL], lambda[efptVDW],
984 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)) != 0,
991 launchPmeGpuSpread(fr->pmedata, box, as_rvec_array(x.unpaddedArrayRef().data()), flags, pmeFlags, wcycle);
994 /* do gridding for pair search */
997 if (graph && bStateChanged)
999 /* Calculate intramolecular shift vectors to make molecules whole */
1000 mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
1004 // - vzero is constant, do we need to pass it?
1005 // - box_diag should be passed directly to nbnxn_put_on_grid
1011 box_diag[XX] = box[XX][XX];
1012 box_diag[YY] = box[YY][YY];
1013 box_diag[ZZ] = box[ZZ][ZZ];
1015 wallcycle_start(wcycle, ewcNS);
1016 if (!DOMAINDECOMP(cr))
1018 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1019 nbnxn_put_on_grid(nbv, box,
1021 nullptr, 0, mdatoms->homenr, -1,
1022 fr->cginfo, x.unpaddedArrayRef(),
1024 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1028 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1029 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd),
1030 fr->cginfo, x.unpaddedArrayRef());
1031 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1034 nbv->setAtomProperties(*mdatoms, fr->cginfo);
1036 wallcycle_stop(wcycle, ewcNS);
1038 /* initialize the GPU nbnxm atom data and bonded data structures */
1041 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1043 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1044 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1045 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1049 /* Now we put all atoms on the grid, we can assign bonded
1050 * interactions to the GPU, where the grid order is
1051 * needed. Also the xq, f and fshift device buffers have
1052 * been reallocated if needed, so the bonded code can
1053 * learn about them. */
1054 // TODO the xq, f, and fshift buffers are now shared
1055 // resources, so they should be maintained by a
1056 // higher-level object than the nb module.
1057 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1059 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1060 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1061 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1063 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1066 // Need to run after the GPU-offload bonded interaction lists
1067 // are set up to be able to determine whether there is bonded work.
1068 setupForceWorkload(ppForceWorkload,
1078 /* do local pair search */
1081 // TODO: fuse this branch with the above bNS block
1082 wallcycle_start_nocount(wcycle, ewcNS);
1083 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1084 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1085 nbv->constructPairlist(Nbnxm::InteractionLocality::Local,
1086 &top->excls, step, nrnb);
1088 nbv->setupGpuShortRangeWork(fr->gpuBonded, Nbnxm::InteractionLocality::Local);
1090 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1091 wallcycle_stop(wcycle, ewcNS);
1093 if (useGpuXBufOps == BufferOpsUseGpu::True)
1095 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1097 // For force buffer ops, we use the below conditon rather than
1098 // useGpuFBufOps to ensure that init is performed even if this
1099 // NS step is also a virial step (on which f buf ops are deactivated).
1100 if (c_enableGpuBufOps && bUseGPU && (GMX_GPU == GMX_GPU_CUDA))
1102 nbv->atomdata_init_add_nbat_f_to_f_gpu();
1107 nbv->setCoordinates(Nbnxm::AtomLocality::Local, false,
1108 x.unpaddedArrayRef(), useGpuXBufOps, pme_gpu_get_device_x(fr->pmedata));
1113 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1115 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1117 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1118 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1119 if (bNS || (useGpuXBufOps == BufferOpsUseGpu::False))
1121 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1122 Nbnxm::AtomLocality::Local);
1124 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1125 // with X buffer ops offloaded to the GPU on all but the search steps
1127 // bonded work not split into separate local and non-local, so with DD
1128 // we can only launch the kernel after non-local coordinates have been received.
1129 if (ppForceWorkload->haveGpuBondedWork && !havePPDomainDecomposition(cr))
1131 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1132 fr->gpuBonded->launchKernel(fr, flags, box);
1133 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1136 /* launch local nonbonded work on GPU */
1137 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1138 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFNo,
1139 step, nrnb, wcycle);
1140 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1141 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1146 // In PME GPU and mixed mode we launch FFT / gather after the
1147 // X copy/transform to allow overlap as well as after the GPU NB
1148 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1149 // the nonbonded kernel.
1150 launchPmeGpuFftAndGather(fr->pmedata, wcycle, useGpuFPmeReduction);
1153 /* Communicate coordinates and sum dipole if necessary +
1154 do non-local pair search */
1155 if (havePPDomainDecomposition(cr))
1159 // TODO: fuse this branch with the above large bNS block
1160 wallcycle_start_nocount(wcycle, ewcNS);
1161 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1162 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1163 nbv->constructPairlist(Nbnxm::InteractionLocality::NonLocal,
1164 &top->excls, step, nrnb);
1166 nbv->setupGpuShortRangeWork(fr->gpuBonded, Nbnxm::InteractionLocality::NonLocal);
1167 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1168 wallcycle_stop(wcycle, ewcNS);
1172 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1174 nbv->setCoordinates(Nbnxm::AtomLocality::NonLocal, false,
1175 x.unpaddedArrayRef(), useGpuXBufOps, pme_gpu_get_device_x(fr->pmedata));
1181 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1183 if (bNS || (useGpuXBufOps == BufferOpsUseGpu::False))
1185 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1186 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1187 Nbnxm::AtomLocality::NonLocal);
1188 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1191 if (ppForceWorkload->haveGpuBondedWork)
1193 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1194 fr->gpuBonded->launchKernel(fr, flags, box);
1195 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1198 /* launch non-local nonbonded tasks on GPU */
1199 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1200 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1201 step, nrnb, wcycle);
1202 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1204 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1210 /* launch D2H copy-back F */
1211 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1212 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1214 bool copyBackNbForce = (useGpuFBufOps == BufferOpsUseGpu::False);
1216 if (havePPDomainDecomposition(cr))
1218 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1219 flags, Nbnxm::AtomLocality::NonLocal, copyBackNbForce);
1221 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1222 flags, Nbnxm::AtomLocality::Local, copyBackNbForce);
1224 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1226 if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1228 fr->gpuBonded->launchEnergyTransfer();
1230 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1233 if (bStateChanged && inputrecNeedMutot(inputrec))
1237 gmx_sumd(2*DIM, mu, cr);
1239 ddBalanceRegionHandler.reopenRegionCpu();
1242 for (i = 0; i < 2; i++)
1244 for (j = 0; j < DIM; j++)
1246 fr->mu_tot[i][j] = mu[i*DIM + j];
1250 if (fr->efep == efepNO)
1252 copy_rvec(fr->mu_tot[0], mu_tot);
1256 for (j = 0; j < DIM; j++)
1259 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1260 lambda[efptCOUL]*fr->mu_tot[1][j];
1264 /* Reset energies */
1265 reset_enerdata(enerd);
1266 clear_rvecs(SHIFTS, fr->fshift);
1268 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1270 wallcycle_start(wcycle, ewcPPDURINGPME);
1271 dd_force_flop_start(cr->dd, nrnb);
1276 wallcycle_start(wcycle, ewcROT);
1277 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, bNS);
1278 wallcycle_stop(wcycle, ewcROT);
1281 /* Start the force cycle counter.
1282 * Note that a different counter is used for dynamic load balancing.
1284 wallcycle_start(wcycle, ewcFORCE);
1286 // set up and clear force outputs
1287 ForceOutputs forceOut = setupForceOutputs(fr, pull_work, *inputrec, force, bDoForces,
1288 ((flags & GMX_FORCE_VIRIAL) != 0), wcycle);
1290 /* We calculate the non-bonded forces, when done on the CPU, here.
1291 * We do this before calling do_force_lowlevel, because in that
1292 * function, the listed forces are calculated before PME, which
1293 * does communication. With this order, non-bonded and listed
1294 * force calculation imbalance can be balanced out by the domain
1295 * decomposition load balancing.
1300 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFYes,
1301 step, nrnb, wcycle);
1304 if (fr->efep != efepNO)
1306 /* Calculate the local and non-local free energy interactions here.
1307 * Happens here on the CPU both with and without GPU.
1309 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::Local,
1310 fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f(), *mdatoms,
1311 inputrec->fepvals, lambda.data(),
1312 enerd, flags, nrnb);
1314 if (havePPDomainDecomposition(cr))
1316 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::NonLocal,
1317 fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f(), *mdatoms,
1318 inputrec->fepvals, lambda.data(),
1319 enerd, flags, nrnb);
1325 if (havePPDomainDecomposition(cr))
1327 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1328 step, nrnb, wcycle);
1331 /* Add all the non-bonded force to the normal force array.
1332 * This can be split into a local and a non-local part when overlapping
1333 * communication with calculation with domain decomposition.
1335 wallcycle_stop(wcycle, ewcFORCE);
1336 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::All, forceOut.f());
1338 wallcycle_start_nocount(wcycle, ewcFORCE);
1340 /* If there are multiple fshift output buffers we need to reduce them */
1341 if (flags & GMX_FORCE_VIRIAL)
1343 /* This is not in a subcounter because it takes a
1344 negligible and constant-sized amount of time */
1345 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->nbat.get(),
1350 /* update QMMMrec, if necessary */
1353 update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1356 /* Compute the bonded and non-bonded energies and optionally forces */
1357 do_force_lowlevel(fr, inputrec, &(top->idef),
1358 cr, ms, nrnb, wcycle, mdatoms,
1359 x, hist, &forceOut, enerd, fcd,
1360 box, lambda.data(), graph, fr->mu_tot,
1362 ddBalanceRegionHandler);
1364 wallcycle_stop(wcycle, ewcFORCE);
1366 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1367 imdSession, pull_work, step, t, wcycle,
1368 fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda.data(),
1369 flags, &forceOut.forceWithVirial(), enerd,
1372 bool useCpuFPmeReduction = thisRankHasDuty(cr, DUTY_PME) && !useGpuFPmeReduction;
1373 bool haveCpuForces = (ppForceWorkload->haveSpecialForces || ppForceWorkload->haveCpuListedForceWork || useCpuFPmeReduction);
1375 // Will store the amount of cycles spent waiting for the GPU that
1376 // will be later used in the DLB accounting.
1377 float cycles_wait_gpu = 0;
1380 /* wait for non-local forces (or calculate in emulation mode) */
1381 if (havePPDomainDecomposition(cr))
1385 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1386 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1387 flags, Nbnxm::AtomLocality::NonLocal,
1388 enerd->grpp.ener[egLJSR].data(),
1389 enerd->grpp.ener[egCOULSR].data(),
1391 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1395 wallcycle_start_nocount(wcycle, ewcFORCE);
1396 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFYes,
1397 step, nrnb, wcycle);
1398 wallcycle_stop(wcycle, ewcFORCE);
1401 if (useGpuFBufOps == BufferOpsUseGpu::True && haveCpuForces)
1403 nbv->launch_copy_f_to_gpu(forceOut.f(), Nbnxm::AtomLocality::NonLocal);
1406 // flag to specify if forces should be accumulated in force buffer
1407 // ops. For non-local part, this just depends on whether CPU forces are present.
1408 bool accumulateForce = (useGpuFBufOps == BufferOpsUseGpu::True) && haveCpuForces;
1409 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::NonLocal,
1410 forceOut.f(), pme_gpu_get_device_f(fr->pmedata),
1411 pme_gpu_get_f_ready_synchronizer(fr->pmedata),
1412 useGpuFBufOps, useGpuFPmeReduction, accumulateForce);
1414 if (useGpuFBufOps == BufferOpsUseGpu::True)
1416 nbv->launch_copy_f_from_gpu(forceOut.f(), Nbnxm::AtomLocality::NonLocal);
1419 if (fr->nbv->emulateGpu() && (flags & GMX_FORCE_VIRIAL))
1421 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->nbat.get(),
1427 if (havePPDomainDecomposition(cr))
1429 /* We are done with the CPU compute.
1430 * We will now communicate the non-local forces.
1431 * If we use a GPU this will overlap with GPU work, so in that case
1432 * we do not close the DD force balancing region here.
1434 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1438 if (useGpuFBufOps == BufferOpsUseGpu::True)
1440 nbv->wait_for_gpu_force_reduction(Nbnxm::AtomLocality::NonLocal);
1442 dd_move_f(cr->dd, force.unpaddedArrayRef(), fr->fshift, wcycle);
1446 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1447 // an alternating wait/reduction scheme.
1448 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPme && bUseGPU && !DOMAINDECOMP(cr) &&
1449 (useGpuFBufOps == BufferOpsUseGpu::False));
1450 if (alternateGpuWait)
1452 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &force, &forceOut.forceWithVirial(), fr->fshift, enerd,
1453 flags, pmeFlags, wcycle);
1456 if (!alternateGpuWait && useGpuPme)
1458 pme_gpu_wait_and_reduce(fr->pmedata, pmeFlags, wcycle, &forceOut.forceWithVirial(), enerd, useGpuFPmeReduction);
1461 /* Wait for local GPU NB outputs on the non-alternating wait path */
1462 if (!alternateGpuWait && bUseGPU)
1464 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1465 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1466 * but even with a step of 0.1 ms the difference is less than 1%
1469 const float gpuWaitApiOverheadMargin = 2e6f; /* cycles */
1471 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1472 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1473 flags, Nbnxm::AtomLocality::Local,
1474 enerd->grpp.ener[egLJSR].data(),
1475 enerd->grpp.ener[egCOULSR].data(),
1477 float cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1479 if (ddBalanceRegionHandler.useBalancingRegion())
1481 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1482 if (bDoForces && cycles_tmp <= gpuWaitApiOverheadMargin)
1484 /* We measured few cycles, it could be that the kernel
1485 * and transfer finished earlier and there was no actual
1486 * wait time, only API call overhead.
1487 * Then the actual time could be anywhere between 0 and
1488 * cycles_wait_est. We will use half of cycles_wait_est.
1490 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1492 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1496 if (fr->nbv->emulateGpu())
1498 // NOTE: emulation kernel is not included in the balancing region,
1499 // but emulation mode does not target performance anyway
1500 wallcycle_start_nocount(wcycle, ewcFORCE);
1501 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local,
1502 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1503 step, nrnb, wcycle);
1504 wallcycle_stop(wcycle, ewcFORCE);
1507 /* Do the nonbonded GPU (or emulation) force buffer reduction
1508 * on the non-alternating path. */
1509 if (bUseOrEmulGPU && !alternateGpuWait)
1512 // TODO: move these steps as early as possible:
1513 // - CPU f H2D should be as soon as all CPU-side forces are done
1514 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
1515 // before the next CPU task that consumes the forces: vsite spread or update)
1517 if (useGpuFBufOps == BufferOpsUseGpu::True && (haveCpuForces || DOMAINDECOMP(cr)))
1519 nbv->launch_copy_f_to_gpu(forceOut.f(), Nbnxm::AtomLocality::Local);
1521 // flag to specify if forces should be accumulated in force
1522 // buffer ops. For local part, this depends on whether CPU
1523 // forces are present, or if DD is active (in which case the
1524 // halo exchange has resulted in contributions from the
1526 bool accumulateForce = (useGpuFBufOps == BufferOpsUseGpu::True) &&
1527 (haveCpuForces || DOMAINDECOMP(cr));
1528 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
1529 forceOut.f(), pme_gpu_get_device_f(fr->pmedata),
1530 pme_gpu_get_f_ready_synchronizer(fr->pmedata),
1531 useGpuFBufOps, useGpuFPmeReduction, accumulateForce);
1533 if (useGpuFBufOps == BufferOpsUseGpu::True)
1535 nbv->launch_copy_f_from_gpu(forceOut.f(), Nbnxm::AtomLocality::Local);
1536 nbv->wait_for_gpu_force_reduction(Nbnxm::AtomLocality::Local);
1540 launchGpuEndOfStepTasks(nbv, fr->gpuBonded, fr->pmedata, enerd,
1547 if (DOMAINDECOMP(cr))
1549 dd_force_flop_stop(cr->dd, nrnb);
1554 /* If we have NoVirSum forces, but we do not calculate the virial,
1555 * we sum fr->f_novirsum=forceOut.f later.
1557 if (vsite && !(fr->haveDirectVirialContributions && !(flags & GMX_FORCE_VIRIAL)))
1559 spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f(), fr->fshift, FALSE, nullptr, nrnb,
1560 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1563 if (flags & GMX_FORCE_VIRIAL)
1565 /* Calculation of the virial must be done after vsites! */
1566 calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f(),
1567 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1571 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1573 /* In case of node-splitting, the PP nodes receive the long-range
1574 * forces, virial and energy from the PME nodes here.
1576 pme_receive_force_ener(cr, &forceOut.forceWithVirial(), enerd, wcycle);
1581 post_process_forces(cr, step, nrnb, wcycle,
1582 top, box, as_rvec_array(x.unpaddedArrayRef().data()), &forceOut,
1583 vir_force, mdatoms, graph, fr, vsite,
1587 if (flags & GMX_FORCE_ENERGY)
1589 /* Sum the potential energy terms from group contributions */
1590 sum_epot(&(enerd->grpp), enerd->term);
1592 if (!EI_TPI(inputrec->eI))
1594 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1598 /* In case we don't have constraints and are using GPUs, the next balancing
1599 * region starts here.
1600 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1601 * virial calculation and COM pulling, is not thus not included in
1602 * the balance timing, which is ok as most tasks do communication.
1604 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);