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);
124 static void sum_forces(rvec f[], gmx::ArrayRef<const gmx::RVec> forceToAdd)
126 const int end = forceToAdd.size();
128 int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
129 #pragma omp parallel for num_threads(nt) schedule(static)
130 for (int i = 0; i < end; i++)
132 rvec_inc(f[i], forceToAdd[i]);
136 static void calc_virial(int start, int homenr, const rvec x[], const rvec f[],
137 tensor vir_part, const t_graph *graph, const matrix box,
138 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
140 /* The short-range virial from surrounding boxes */
141 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
142 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
144 /* Calculate partial virial, for local atoms only, based on short range.
145 * Total virial is computed in global_stat, called from do_md
147 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
148 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
152 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
156 static void pull_potential_wrapper(const t_commrec *cr,
157 const t_inputrec *ir,
158 const matrix box, gmx::ArrayRef<const gmx::RVec> x,
159 gmx::ForceWithVirial *force,
160 const t_mdatoms *mdatoms,
161 gmx_enerdata_t *enerd,
165 gmx_wallcycle_t wcycle)
170 /* Calculate the center of mass forces, this requires communication,
171 * which is why pull_potential is called close to other communication.
173 wallcycle_start(wcycle, ewcPULLPOT);
174 set_pbc(&pbc, ir->ePBC, box);
176 enerd->term[F_COM_PULL] +=
177 pull_potential(pull_work, mdatoms, &pbc,
178 cr, t, lambda[efptRESTRAINT], as_rvec_array(x.data()), force, &dvdl);
179 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
180 wallcycle_stop(wcycle, ewcPULLPOT);
183 static void pme_receive_force_ener(const t_commrec *cr,
184 gmx::ForceWithVirial *forceWithVirial,
185 gmx_enerdata_t *enerd,
186 gmx_wallcycle_t wcycle)
188 real e_q, e_lj, dvdl_q, dvdl_lj;
189 float cycles_ppdpme, cycles_seppme;
191 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
192 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
194 /* In case of node-splitting, the PP nodes receive the long-range
195 * forces, virial and energy from the PME nodes here.
197 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
200 gmx_pme_receive_f(cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
202 enerd->term[F_COUL_RECIP] += e_q;
203 enerd->term[F_LJ_RECIP] += e_lj;
204 enerd->dvdl_lin[efptCOUL] += dvdl_q;
205 enerd->dvdl_lin[efptVDW] += dvdl_lj;
209 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
211 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
214 static void print_large_forces(FILE *fp,
222 real force2Tolerance = gmx::square(forceTolerance);
223 gmx::index numNonFinite = 0;
224 for (int i = 0; i < md->homenr; i++)
226 real force2 = norm2(f[i]);
227 bool nonFinite = !std::isfinite(force2);
228 if (force2 >= force2Tolerance || nonFinite)
230 fprintf(fp, "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
232 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
239 if (numNonFinite > 0)
241 /* Note that with MPI this fatal call on one rank might interrupt
242 * the printing on other ranks. But we can only avoid that with
243 * an expensive MPI barrier that we would need at each step.
245 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
249 static void post_process_forces(const t_commrec *cr,
252 gmx_wallcycle_t wcycle,
253 const gmx_localtop_t *top,
257 gmx::ForceWithVirial *forceWithVirial,
259 const t_mdatoms *mdatoms,
260 const t_graph *graph,
261 const t_forcerec *fr,
262 const gmx_vsite_t *vsite,
265 if (fr->haveDirectVirialContributions)
267 rvec *fDirectVir = as_rvec_array(forceWithVirial->force_.data());
271 /* Spread the mesh force on virtual sites to the other particles...
272 * This is parallellized. MPI communication is performed
273 * if the constructing atoms aren't local.
275 matrix virial = { { 0 } };
276 spread_vsite_f(vsite, x, fDirectVir, nullptr,
277 (flags & GMX_FORCE_VIRIAL) != 0, virial,
279 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
280 forceWithVirial->addVirialContribution(virial);
283 if (flags & GMX_FORCE_VIRIAL)
285 /* Now add the forces, this is local */
286 sum_forces(f, forceWithVirial->force_);
288 /* Add the direct virial contributions */
289 GMX_ASSERT(forceWithVirial->computeVirial_, "forceWithVirial should request virial computation when we request the virial");
290 m_add(vir_force, forceWithVirial->getVirial(), vir_force);
294 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
299 if (fr->print_force >= 0)
301 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
305 static void do_nb_verlet(t_forcerec *fr,
306 const interaction_const_t *ic,
307 gmx_enerdata_t *enerd,
309 const Nbnxm::InteractionLocality ilocality,
313 gmx_wallcycle_t wcycle)
315 if (!(flags & GMX_FORCE_NONBONDED))
317 /* skip non-bonded calculation */
321 nonbonded_verlet_t *nbv = fr->nbv.get();
323 /* GPU kernel launch overhead is already timed separately */
324 if (fr->cutoff_scheme != ecutsVERLET)
326 gmx_incons("Invalid cut-off scheme passed!");
331 /* When dynamic pair-list pruning is requested, we need to prune
332 * at nstlistPrune steps.
334 if (nbv->isDynamicPruningStepCpu(step))
336 /* Prune the pair-list beyond fr->ic->rlistPrune using
337 * the current coordinates of the atoms.
339 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
340 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
341 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
344 wallcycle_sub_start(wcycle, ewcsNONBONDED);
347 nbv->dispatchNonbondedKernel(ilocality, *ic, flags, clearF, fr, enerd, nrnb);
351 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
355 static inline void clear_rvecs_omp(int n, rvec v[])
357 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, n);
359 /* Note that we would like to avoid this conditional by putting it
360 * into the omp pragma instead, but then we still take the full
361 * omp parallel for overhead (at least with gcc5).
365 for (int i = 0; i < n; i++)
372 #pragma omp parallel for num_threads(nth) schedule(static)
373 for (int i = 0; i < n; i++)
380 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
382 * \param groupOptions Group options, containing T-coupling options
384 static real averageKineticEnergyEstimate(const t_grpopts &groupOptions)
386 real nrdfCoupled = 0;
387 real nrdfUncoupled = 0;
388 real kineticEnergy = 0;
389 for (int g = 0; g < groupOptions.ngtc; g++)
391 if (groupOptions.tau_t[g] >= 0)
393 nrdfCoupled += groupOptions.nrdf[g];
394 kineticEnergy += groupOptions.nrdf[g]*0.5*groupOptions.ref_t[g]*BOLTZ;
398 nrdfUncoupled += groupOptions.nrdf[g];
402 /* This conditional with > also catches nrdf=0 */
403 if (nrdfCoupled > nrdfUncoupled)
405 return kineticEnergy*(nrdfCoupled + nrdfUncoupled)/nrdfCoupled;
413 /*! \brief This routine checks that the potential energy is finite.
415 * Always checks that the potential energy is finite. If step equals
416 * inputrec.init_step also checks that the magnitude of the potential energy
417 * is reasonable. Terminates with a fatal error when a check fails.
418 * Note that passing this check does not guarantee finite forces,
419 * since those use slightly different arithmetics. But in most cases
420 * there is just a narrow coordinate range where forces are not finite
421 * and energies are finite.
423 * \param[in] step The step number, used for checking and printing
424 * \param[in] enerd The energy data; the non-bonded group energies need to be added to enerd.term[F_EPOT] before calling this routine
425 * \param[in] inputrec The input record
427 static void checkPotentialEnergyValidity(int64_t step,
428 const gmx_enerdata_t &enerd,
429 const t_inputrec &inputrec)
431 /* Threshold valid for comparing absolute potential energy against
432 * the kinetic energy. Normally one should not consider absolute
433 * potential energy values, but with a factor of one million
434 * we should never get false positives.
436 constexpr real c_thresholdFactor = 1e6;
438 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
439 real averageKineticEnergy = 0;
440 /* We only check for large potential energy at the initial step,
441 * because that is by far the most likely step for this too occur
442 * and because computing the average kinetic energy is not free.
443 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
444 * before they become NaN.
446 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
448 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
451 if (energyIsNotFinite || (averageKineticEnergy > 0 &&
452 enerd.term[F_EPOT] > c_thresholdFactor*averageKineticEnergy))
454 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.",
457 energyIsNotFinite ? "not finite" : "extremely high",
459 enerd.term[F_COUL_SR],
460 energyIsNotFinite ? "non-finite" : "very high",
461 energyIsNotFinite ? " or Nan" : "");
465 /*! \brief Return true if there are special forces computed this step.
467 * The conditionals exactly correspond to those in computeSpecialForces().
470 haveSpecialForces(const t_inputrec *inputrec,
471 ForceProviders *forceProviders,
472 const pull_t *pull_work,
476 const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
479 ((computeForces && forceProviders->hasForceProvider()) || // forceProviders
480 (inputrec->bPull && pull_have_potential(pull_work)) || // pull
481 inputrec->bRot || // enforced rotation
482 (ed != nullptr) || // flooding
483 (inputrec->bIMD && computeForces)); // IMD
486 /*! \brief Compute forces and/or energies for special algorithms
488 * The intention is to collect all calls to algorithms that compute
489 * forces on local atoms only and that do not contribute to the local
490 * virial sum (but add their virial contribution separately).
491 * Eventually these should likely all become ForceProviders.
492 * Within this function the intention is to have algorithms that do
493 * global communication at the end, so global barriers within the MD loop
494 * are as close together as possible.
496 * \param[in] fplog The log file
497 * \param[in] cr The communication record
498 * \param[in] inputrec The input record
499 * \param[in] awh The Awh module (nullptr if none in use).
500 * \param[in] enforcedRotation Enforced rotation module.
501 * \param[in] imdSession The IMD session
502 * \param[in] pull_work The pull work structure.
503 * \param[in] step The current MD step
504 * \param[in] t The current time
505 * \param[in,out] wcycle Wallcycle accounting struct
506 * \param[in,out] forceProviders Pointer to a list of force providers
507 * \param[in] box The unit cell
508 * \param[in] x The coordinates
509 * \param[in] mdatoms Per atom properties
510 * \param[in] lambda Array of free-energy lambda values
511 * \param[in] forceFlags Flags that tell whether we should compute forces/energies/virial
512 * \param[in,out] forceWithVirial Force and virial buffers
513 * \param[in,out] enerd Energy buffer
514 * \param[in,out] ed Essential dynamics pointer
515 * \param[in] bNS Tells if we did neighbor searching this step, used for ED sampling
517 * \todo Remove bNS, which is used incorrectly.
518 * \todo Convert all other algorithms called here to ForceProviders.
521 computeSpecialForces(FILE *fplog,
523 const t_inputrec *inputrec,
525 gmx_enfrot *enforcedRotation,
526 gmx::ImdSession *imdSession,
530 gmx_wallcycle_t wcycle,
531 ForceProviders *forceProviders,
533 gmx::ArrayRef<const gmx::RVec> x,
534 const t_mdatoms *mdatoms,
537 gmx::ForceWithVirial *forceWithVirial,
538 gmx_enerdata_t *enerd,
542 const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
544 /* NOTE: Currently all ForceProviders only provide forces.
545 * When they also provide energies, remove this conditional.
549 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
550 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
552 /* Collect forces from modules */
553 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
556 if (inputrec->bPull && pull_have_potential(pull_work))
558 pull_potential_wrapper(cr, inputrec, box, x,
560 mdatoms, enerd, pull_work, lambda, t,
565 enerd->term[F_COM_PULL] +=
566 awh->applyBiasForcesAndUpdateBias(inputrec->ePBC, *mdatoms, box,
568 t, step, wcycle, fplog);
572 rvec *f = as_rvec_array(forceWithVirial->force_.data());
574 /* Add the forces from enforced rotation potentials (if any) */
577 wallcycle_start(wcycle, ewcROTadd);
578 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
579 wallcycle_stop(wcycle, ewcROTadd);
584 /* Note that since init_edsam() is called after the initialization
585 * of forcerec, edsam doesn't request the noVirSum force buffer.
586 * Thus if no other algorithm (e.g. PME) requires it, the forces
587 * here will contribute to the virial.
589 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, bNS);
592 /* Add forces from interactive molecular dynamics (IMD), if any */
593 if (inputrec->bIMD && computeForces)
595 imdSession->applyForces(f);
599 /*! \brief Launch the prepare_step and spread stages of PME GPU.
601 * \param[in] pmedata The PME structure
602 * \param[in] box The box matrix
603 * \param[in] x Coordinate array
604 * \param[in] flags Force flags
605 * \param[in] pmeFlags PME flags
606 * \param[in] wcycle The wallcycle structure
608 static inline void launchPmeGpuSpread(gmx_pme_t *pmedata,
613 gmx_wallcycle_t wcycle)
615 pme_gpu_prepare_computation(pmedata, (flags & GMX_FORCE_DYNAMICBOX) != 0, box, wcycle, pmeFlags);
616 pme_gpu_launch_spread(pmedata, x, wcycle);
619 /*! \brief Launch the FFT and gather stages of PME GPU
621 * This function only implements setting the output forces (no accumulation).
623 * \param[in] pmedata The PME structure
624 * \param[in] wcycle The wallcycle structure
626 static void launchPmeGpuFftAndGather(gmx_pme_t *pmedata,
627 gmx_wallcycle_t wcycle)
629 pme_gpu_launch_complex_transforms(pmedata, wcycle);
630 pme_gpu_launch_gather(pmedata, wcycle, PmeForceOutputHandling::Set);
634 * Polling wait for either of the PME or nonbonded GPU tasks.
636 * Instead of a static order in waiting for GPU tasks, this function
637 * polls checking which of the two tasks completes first, and does the
638 * associated force buffer reduction overlapped with the other task.
639 * By doing that, unlike static scheduling order, it can always overlap
640 * one of the reductions, regardless of the GPU task completion order.
642 * \param[in] nbv Nonbonded verlet structure
643 * \param[in,out] pmedata PME module data
644 * \param[in,out] force Force array to reduce task outputs into.
645 * \param[in,out] forceWithVirial Force and virial buffers
646 * \param[in,out] fshift Shift force output vector results are reduced into
647 * \param[in,out] enerd Energy data structure results are reduced into
648 * \param[in] flags Force flags
649 * \param[in] pmeFlags PME flags
650 * \param[in] haveOtherWork Tells whether there is other work than non-bonded in the stream(s)
651 * \param[in] wcycle The wallcycle structure
653 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t *nbv,
655 gmx::ArrayRefWithPadding<gmx::RVec> *force,
656 gmx::ForceWithVirial *forceWithVirial,
658 gmx_enerdata_t *enerd,
662 gmx_wallcycle_t wcycle)
664 bool isPmeGpuDone = false;
665 bool isNbGpuDone = false;
668 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
670 while (!isPmeGpuDone || !isNbGpuDone)
674 GpuTaskCompletion completionType = (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
675 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, pmeFlags, wcycle, forceWithVirial, enerd, completionType);
680 GpuTaskCompletion completionType = (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
681 wallcycle_start_nocount(wcycle, ewcWAIT_GPU_NB_L);
682 isNbGpuDone = Nbnxm::gpu_try_finish_task(nbv->gpu_nbv,
684 Nbnxm::AtomLocality::Local,
686 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
687 fshift, completionType);
688 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
689 // To get the call count right, when the task finished we
690 // issue a start/stop.
691 // TODO: move the ewcWAIT_GPU_NB_L cycle counting into nbnxn_gpu_try_finish_task()
692 // and ewcNB_XF_BUF_OPS counting into nbnxn_atomdata_add_nbat_f_to_f().
695 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
696 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
698 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
699 as_rvec_array(force->unpaddedArrayRef().data()), wcycle);
705 /*! \brief Hack structure with force ouput buffers for do_force for the home atoms for this domain */
709 ForceOutputs(rvec *f, gmx::ForceWithVirial const forceWithVirial) :
711 forceWithVirial(forceWithVirial) {}
713 //! Force output buffer used by legacy modules (without SIMD padding)
715 //! Force with direct virial contribution (if there are any; without SIMD padding)
716 gmx::ForceWithVirial forceWithVirial;
719 /*! \brief Set up the different force buffers; also does clearing.
721 * \param[in] fr force record pointer
722 * \param[in] pull_work The pull work object.
723 * \param[in] inputrec input record
724 * \param[in] force force array
725 * \param[in] bDoForces True if force are computed this step
726 * \param[in] doVirial True if virial is computed this step
727 * \param[out] wcycle wallcycle recording structure
729 * \returns Cleared force output structure
732 setupForceOutputs(const t_forcerec *fr,
734 const t_inputrec &inputrec,
735 gmx::ArrayRefWithPadding<gmx::RVec> force,
736 const bool bDoForces,
738 gmx_wallcycle_t wcycle)
740 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
742 /* Temporary solution until all routines take PaddedRVecVector */
743 rvec *const f = as_rvec_array(force.unpaddedArrayRef().data());
746 /* Clear the short- and long-range forces */
747 clear_rvecs_omp(fr->natoms_force_constr, f);
750 /* If we need to compute the virial, we might need a separate
751 * force buffer for algorithms for which the virial is calculated
752 * directly, such as PME. Otherwise, forceWithVirial uses the
753 * the same force (f in legacy calls) buffer as other algorithms.
755 const bool useSeparateForceWithVirialBuffer = (bDoForces && (doVirial && fr->haveDirectVirialContributions));
758 /* forceWithVirial uses the local atom range only */
759 gmx::ForceWithVirial forceWithVirial (useSeparateForceWithVirialBuffer ?
760 *fr->forceBufferForDirectVirialContributions : force.unpaddedArrayRef(),
763 if (useSeparateForceWithVirialBuffer)
765 /* TODO: update comment
766 * We only compute forces on local atoms. Note that vsites can
767 * spread to non-local atoms, but that part of the buffer is
768 * cleared separately in the vsite spreading code.
770 clear_rvecs_omp(forceWithVirial.force_.size(), as_rvec_array(forceWithVirial.force_.data()));
773 if (inputrec.bPull && pull_have_constraint(pull_work))
775 clear_pull_forces(pull_work);
778 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
780 return ForceOutputs(f, forceWithVirial);
784 /*! \brief Set up flags that indicate what type of work is there to compute.
786 * Currently we only update it at search steps,
787 * but some properties may change more frequently (e.g. virial/non-virial step),
788 * so when including those either the frequency of update (per-step) or the scope
789 * of a flag will change (i.e. a set of flags for nstlist steps).
793 setupForceWorkload(gmx::PpForceWorkload *forceWork,
794 const t_inputrec *inputrec,
795 const t_forcerec *fr,
796 const pull_t *pull_work,
803 forceWork->haveSpecialForces = haveSpecialForces(inputrec, fr->forceProviders, pull_work, forceFlags, ed);
804 forceWork->haveCpuBondedWork = haveCpuBondeds(*fr);
805 forceWork->haveGpuBondedWork = ((fr->gpuBonded != nullptr) && fr->gpuBonded->haveInteractions());
806 forceWork->haveRestraintsWork = havePositionRestraints(idef, *fcd);
807 forceWork->haveCpuListedForceWork = haveCpuListedForces(*fr, idef, *fcd);
810 void do_force(FILE *fplog,
812 const gmx_multisim_t *ms,
813 const t_inputrec *inputrec,
815 gmx_enfrot *enforcedRotation,
816 gmx::ImdSession *imdSession,
820 gmx_wallcycle_t wcycle,
821 const gmx_localtop_t *top,
823 gmx::ArrayRefWithPadding<gmx::RVec> x, //NOLINT(performance-unnecessary-value-param)
825 gmx::ArrayRefWithPadding<gmx::RVec> force, //NOLINT(performance-unnecessary-value-param)
827 const t_mdatoms *mdatoms,
828 gmx_enerdata_t *enerd,
830 gmx::ArrayRef<real> lambda,
833 gmx::PpForceWorkload *ppForceWorkload,
834 const gmx_vsite_t *vsite,
839 const DDBalanceRegionHandler &ddBalanceRegionHandler)
843 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
844 gmx_bool bDoForces, bUseGPU, bUseOrEmulGPU;
845 nonbonded_verlet_t *nbv = fr->nbv.get();
846 interaction_const_t *ic = fr->ic;
848 /* modify force flag if not doing nonbonded */
851 flags &= ~GMX_FORCE_NONBONDED;
853 bStateChanged = ((flags & GMX_FORCE_STATECHANGED) != 0);
854 bNS = ((flags & GMX_FORCE_NS) != 0);
855 bFillGrid = (bNS && bStateChanged);
856 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
857 bDoForces = ((flags & GMX_FORCE_FORCES) != 0);
858 bUseGPU = fr->nbv->useGpu();
859 bUseOrEmulGPU = bUseGPU || fr->nbv->emulateGpu();
861 const auto pmeRunMode = fr->pmedata ? pme_run_mode(fr->pmedata) : PmeRunMode::CPU;
862 // TODO slim this conditional down - inputrec and duty checks should mean the same in proper code!
863 const bool useGpuPme = EEL_PME(fr->ic->eeltype) && thisRankHasDuty(cr, DUTY_PME) &&
864 ((pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed));
865 const int pmeFlags = GMX_PME_SPREAD | GMX_PME_SOLVE |
866 ((flags & GMX_FORCE_VIRIAL) ? GMX_PME_CALC_ENER_VIR : 0) |
867 ((flags & GMX_FORCE_ENERGY) ? GMX_PME_CALC_ENER_VIR : 0) |
868 ((flags & GMX_FORCE_FORCES) ? GMX_PME_CALC_F : 0);
870 /* At a search step we need to start the first balancing region
871 * somewhere early inside the step after communication during domain
872 * decomposition (and not during the previous step as usual).
876 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
880 const int homenr = mdatoms->homenr;
882 clear_mat(vir_force);
886 update_forcerec(fr, box);
888 if (inputrecNeedMutot(inputrec))
890 /* Calculate total (local) dipole moment in a temporary common array.
891 * This makes it possible to sum them over nodes faster.
893 calc_mu(start, homenr,
894 x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
899 if (fr->ePBC != epbcNONE)
901 /* Compute shift vectors every step,
902 * because of pressure coupling or box deformation!
904 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
906 calc_shifts(box, fr->shift_vec);
911 put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr), gmx_omp_nthreads_get(emntDefault));
912 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
914 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
916 unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
920 nbnxn_atomdata_copy_shiftvec((flags & GMX_FORCE_DYNAMICBOX) != 0,
921 fr->shift_vec, nbv->nbat.get());
924 if (!thisRankHasDuty(cr, DUTY_PME))
926 /* Send particle coordinates to the pme nodes.
927 * Since this is only implemented for domain decomposition
928 * and domain decomposition does not use the graph,
929 * we do not need to worry about shifting.
931 gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
932 lambda[efptCOUL], lambda[efptVDW],
933 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)) != 0,
940 launchPmeGpuSpread(fr->pmedata, box, as_rvec_array(x.unpaddedArrayRef().data()), flags, pmeFlags, wcycle);
943 /* do gridding for pair search */
946 if (graph && bStateChanged)
948 /* Calculate intramolecular shift vectors to make molecules whole */
949 mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
953 // - vzero is constant, do we need to pass it?
954 // - box_diag should be passed directly to nbnxn_put_on_grid
960 box_diag[XX] = box[XX][XX];
961 box_diag[YY] = box[YY][YY];
962 box_diag[ZZ] = box[ZZ][ZZ];
964 wallcycle_start(wcycle, ewcNS);
965 if (!DOMAINDECOMP(cr))
967 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
968 nbnxn_put_on_grid(nbv, box,
970 nullptr, 0, mdatoms->homenr, -1,
971 fr->cginfo, x.unpaddedArrayRef(),
973 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
977 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
978 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd),
979 fr->cginfo, x.unpaddedArrayRef());
980 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
983 nbv->setAtomProperties(*mdatoms, *fr->cginfo);
985 wallcycle_stop(wcycle, ewcNS);
987 /* initialize the GPU nbnxm atom data and bonded data structures */
990 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
992 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
993 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
994 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
998 /* Now we put all atoms on the grid, we can assign bonded
999 * interactions to the GPU, where the grid order is
1000 * needed. Also the xq, f and fshift device buffers have
1001 * been reallocated if needed, so the bonded code can
1002 * learn about them. */
1003 // TODO the xq, f, and fshift buffers are now shared
1004 // resources, so they should be maintained by a
1005 // higher-level object than the nb module.
1006 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1008 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1009 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1010 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1012 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1015 // Need to run after the GPU-offload bonded interaction lists
1016 // are set up to be able to determine whether there is bonded work.
1017 setupForceWorkload(ppForceWorkload,
1027 /* do local pair search */
1030 // TODO: fuse this branch with the above bNS block
1031 wallcycle_start_nocount(wcycle, ewcNS);
1032 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1033 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1034 nbv->constructPairlist(Nbnxm::InteractionLocality::Local,
1035 &top->excls, step, nrnb);
1036 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1037 wallcycle_stop(wcycle, ewcNS);
1041 nbv->setCoordinates(Nbnxm::AtomLocality::Local, false,
1042 x.unpaddedArrayRef(), wcycle);
1047 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1049 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1051 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1052 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1053 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1054 Nbnxm::AtomLocality::Local,
1055 ppForceWorkload->haveGpuBondedWork);
1056 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1058 // bonded work not split into separate local and non-local, so with DD
1059 // we can only launch the kernel after non-local coordinates have been received.
1060 if (ppForceWorkload->haveGpuBondedWork && !havePPDomainDecomposition(cr))
1062 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1063 fr->gpuBonded->launchKernels(fr, flags, box);
1064 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1067 /* launch local nonbonded work on GPU */
1068 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1069 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFNo,
1070 step, nrnb, wcycle);
1071 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1072 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1077 // In PME GPU and mixed mode we launch FFT / gather after the
1078 // X copy/transform to allow overlap as well as after the GPU NB
1079 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1080 // the nonbonded kernel.
1081 launchPmeGpuFftAndGather(fr->pmedata, wcycle);
1084 /* Communicate coordinates and sum dipole if necessary +
1085 do non-local pair search */
1086 if (havePPDomainDecomposition(cr))
1090 // TODO: fuse this branch with the above large bNS block
1091 wallcycle_start_nocount(wcycle, ewcNS);
1092 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1093 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1094 nbv->constructPairlist(Nbnxm::InteractionLocality::NonLocal,
1095 &top->excls, step, nrnb);
1096 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1097 wallcycle_stop(wcycle, ewcNS);
1101 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1103 nbv->setCoordinates(Nbnxm::AtomLocality::NonLocal, false,
1104 x.unpaddedArrayRef(), wcycle);
1109 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1111 /* launch non-local nonbonded tasks on GPU */
1112 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1113 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1114 Nbnxm::AtomLocality::NonLocal,
1115 ppForceWorkload->haveGpuBondedWork);
1116 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1118 if (ppForceWorkload->haveGpuBondedWork)
1120 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1121 fr->gpuBonded->launchKernels(fr, flags, box);
1122 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1125 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1126 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1127 step, nrnb, wcycle);
1128 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1130 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1136 /* launch D2H copy-back F */
1137 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1138 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1139 if (havePPDomainDecomposition(cr))
1141 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1142 flags, Nbnxm::AtomLocality::NonLocal, ppForceWorkload->haveGpuBondedWork);
1144 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1145 flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork);
1146 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1148 if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1150 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1151 fr->gpuBonded->launchEnergyTransfer();
1152 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1154 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1157 if (bStateChanged && inputrecNeedMutot(inputrec))
1161 gmx_sumd(2*DIM, mu, cr);
1163 ddBalanceRegionHandler.reopenRegionCpu();
1166 for (i = 0; i < 2; i++)
1168 for (j = 0; j < DIM; j++)
1170 fr->mu_tot[i][j] = mu[i*DIM + j];
1174 if (fr->efep == efepNO)
1176 copy_rvec(fr->mu_tot[0], mu_tot);
1180 for (j = 0; j < DIM; j++)
1183 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1184 lambda[efptCOUL]*fr->mu_tot[1][j];
1188 /* Reset energies */
1189 reset_enerdata(enerd);
1190 clear_rvecs(SHIFTS, fr->fshift);
1192 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1194 wallcycle_start(wcycle, ewcPPDURINGPME);
1195 dd_force_flop_start(cr->dd, nrnb);
1200 wallcycle_start(wcycle, ewcROT);
1201 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, bNS);
1202 wallcycle_stop(wcycle, ewcROT);
1205 /* Start the force cycle counter.
1206 * Note that a different counter is used for dynamic load balancing.
1208 wallcycle_start(wcycle, ewcFORCE);
1210 // set up and clear force outputs
1211 struct ForceOutputs forceOut = setupForceOutputs(fr, pull_work, *inputrec, force, bDoForces,
1212 ((flags & GMX_FORCE_VIRIAL) != 0), wcycle);
1214 /* We calculate the non-bonded forces, when done on the CPU, here.
1215 * We do this before calling do_force_lowlevel, because in that
1216 * function, the listed forces are calculated before PME, which
1217 * does communication. With this order, non-bonded and listed
1218 * force calculation imbalance can be balanced out by the domain
1219 * decomposition load balancing.
1224 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFYes,
1225 step, nrnb, wcycle);
1228 if (fr->efep != efepNO)
1230 /* Calculate the local and non-local free energy interactions here.
1231 * Happens here on the CPU both with and without GPU.
1233 wallcycle_sub_start(wcycle, ewcsNONBONDED);
1234 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::Local,
1235 fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1236 inputrec->fepvals, lambda.data(),
1237 enerd, flags, nrnb);
1239 if (havePPDomainDecomposition(cr))
1241 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::NonLocal,
1242 fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1243 inputrec->fepvals, lambda.data(),
1244 enerd, flags, nrnb);
1246 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
1251 if (havePPDomainDecomposition(cr))
1253 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1254 step, nrnb, wcycle);
1257 /* Add all the non-bonded force to the normal force array.
1258 * This can be split into a local and a non-local part when overlapping
1259 * communication with calculation with domain decomposition.
1261 wallcycle_stop(wcycle, ewcFORCE);
1263 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::All, forceOut.f, wcycle);
1265 wallcycle_start_nocount(wcycle, ewcFORCE);
1267 /* If there are multiple fshift output buffers we need to reduce them */
1268 if (flags & GMX_FORCE_VIRIAL)
1270 /* This is not in a subcounter because it takes a
1271 negligible and constant-sized amount of time */
1272 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->nbat.get(),
1277 /* update QMMMrec, if necessary */
1280 update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1283 /* Compute the bonded and non-bonded energies and optionally forces */
1284 do_force_lowlevel(fr, inputrec, &(top->idef),
1285 cr, ms, nrnb, wcycle, mdatoms,
1286 as_rvec_array(x.unpaddedArrayRef().data()), hist, forceOut.f, &forceOut.forceWithVirial, enerd, fcd,
1287 box, lambda.data(), graph, &(top->excls), fr->mu_tot,
1289 ddBalanceRegionHandler);
1291 wallcycle_stop(wcycle, ewcFORCE);
1293 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1294 imdSession, pull_work, step, t, wcycle,
1295 fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda.data(),
1296 flags, &forceOut.forceWithVirial, enerd,
1299 // Will store the amount of cycles spent waiting for the GPU that
1300 // will be later used in the DLB accounting.
1301 float cycles_wait_gpu = 0;
1304 /* wait for non-local forces (or calculate in emulation mode) */
1305 if (havePPDomainDecomposition(cr))
1309 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1310 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1311 flags, Nbnxm::AtomLocality::NonLocal,
1312 ppForceWorkload->haveGpuBondedWork,
1313 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1315 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1319 wallcycle_start_nocount(wcycle, ewcFORCE);
1320 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFYes,
1321 step, nrnb, wcycle);
1322 wallcycle_stop(wcycle, ewcFORCE);
1325 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::NonLocal,
1326 forceOut.f, wcycle);
1330 if (havePPDomainDecomposition(cr))
1332 /* We are done with the CPU compute.
1333 * We will now communicate the non-local forces.
1334 * If we use a GPU this will overlap with GPU work, so in that case
1335 * we do not close the DD force balancing region here.
1337 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1341 dd_move_f(cr->dd, force.unpaddedArrayRef(), fr->fshift, wcycle);
1345 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1346 // an alternating wait/reduction scheme.
1347 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPme && bUseGPU && !DOMAINDECOMP(cr));
1348 if (alternateGpuWait)
1350 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &force, &forceOut.forceWithVirial, fr->fshift, enerd,
1351 flags, pmeFlags, ppForceWorkload->haveGpuBondedWork, wcycle);
1354 if (!alternateGpuWait && useGpuPme)
1356 pme_gpu_wait_and_reduce(fr->pmedata, pmeFlags, wcycle, &forceOut.forceWithVirial, enerd);
1359 /* Wait for local GPU NB outputs on the non-alternating wait path */
1360 if (!alternateGpuWait && bUseGPU)
1362 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1363 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1364 * but even with a step of 0.1 ms the difference is less than 1%
1367 const float gpuWaitApiOverheadMargin = 2e6f; /* cycles */
1369 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1370 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1371 flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork,
1372 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1374 float cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1376 if (ddBalanceRegionHandler.useBalancingRegion())
1378 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1379 if (bDoForces && cycles_tmp <= gpuWaitApiOverheadMargin)
1381 /* We measured few cycles, it could be that the kernel
1382 * and transfer finished earlier and there was no actual
1383 * wait time, only API call overhead.
1384 * Then the actual time could be anywhere between 0 and
1385 * cycles_wait_est. We will use half of cycles_wait_est.
1387 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1389 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1393 if (fr->nbv->emulateGpu())
1395 // NOTE: emulation kernel is not included in the balancing region,
1396 // but emulation mode does not target performance anyway
1397 wallcycle_start_nocount(wcycle, ewcFORCE);
1398 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local,
1399 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1400 step, nrnb, wcycle);
1401 wallcycle_stop(wcycle, ewcFORCE);
1406 pme_gpu_reinit_computation(fr->pmedata, wcycle);
1411 /* now clear the GPU outputs while we finish the step on the CPU */
1412 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1413 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1414 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, flags);
1416 if (nbv->isDynamicPruningStepGpu(step))
1418 nbv->dispatchPruneKernelGpu(step);
1420 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1421 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1424 if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1426 wallcycle_start(wcycle, ewcWAIT_GPU_BONDED);
1427 // in principle this should be included in the DD balancing region,
1428 // but generally it is infrequent so we'll omit it for the sake of
1430 fr->gpuBonded->accumulateEnergyTerms(enerd);
1431 wallcycle_stop(wcycle, ewcWAIT_GPU_BONDED);
1433 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1434 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1435 fr->gpuBonded->clearEnergies();
1436 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1437 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1440 /* Do the nonbonded GPU (or emulation) force buffer reduction
1441 * on the non-alternating path. */
1442 if (bUseOrEmulGPU && !alternateGpuWait)
1444 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
1445 forceOut.f, wcycle);
1447 if (DOMAINDECOMP(cr))
1449 dd_force_flop_stop(cr->dd, nrnb);
1454 /* If we have NoVirSum forces, but we do not calculate the virial,
1455 * we sum fr->f_novirsum=forceOut.f later.
1457 if (vsite && !(fr->haveDirectVirialContributions && !(flags & GMX_FORCE_VIRIAL)))
1459 spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, fr->fshift, FALSE, nullptr, nrnb,
1460 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1463 if (flags & GMX_FORCE_VIRIAL)
1465 /* Calculation of the virial must be done after vsites! */
1466 calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f,
1467 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1471 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1473 /* In case of node-splitting, the PP nodes receive the long-range
1474 * forces, virial and energy from the PME nodes here.
1476 pme_receive_force_ener(cr, &forceOut.forceWithVirial, enerd, wcycle);
1481 post_process_forces(cr, step, nrnb, wcycle,
1482 top, box, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, &forceOut.forceWithVirial,
1483 vir_force, mdatoms, graph, fr, vsite,
1487 if (flags & GMX_FORCE_ENERGY)
1489 /* Sum the potential energy terms from group contributions */
1490 sum_epot(&(enerd->grpp), enerd->term);
1492 if (!EI_TPI(inputrec->eI))
1494 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1498 /* In case we don't have constraints and are using GPUs, the next balancing
1499 * region starts here.
1500 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1501 * virial calculation and COM pulling, is not thus not included in
1502 * the balance timing, which is ok as most tasks do communication.
1504 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);