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/force.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/gmx_omp_nthreads.h"
80 #include "gromacs/mdlib/ppforceworkload.h"
81 #include "gromacs/mdlib/qmmm.h"
82 #include "gromacs/mdlib/update.h"
83 #include "gromacs/mdtypes/commrec.h"
84 #include "gromacs/mdtypes/enerdata.h"
85 #include "gromacs/mdtypes/forceoutput.h"
86 #include "gromacs/mdtypes/iforceprovider.h"
87 #include "gromacs/mdtypes/inputrec.h"
88 #include "gromacs/mdtypes/md_enums.h"
89 #include "gromacs/mdtypes/state.h"
90 #include "gromacs/nbnxm/atomdata.h"
91 #include "gromacs/nbnxm/gpu_data_mgmt.h"
92 #include "gromacs/nbnxm/nbnxm.h"
93 #include "gromacs/pbcutil/ishift.h"
94 #include "gromacs/pbcutil/mshift.h"
95 #include "gromacs/pbcutil/pbc.h"
96 #include "gromacs/pulling/pull.h"
97 #include "gromacs/pulling/pull_rotation.h"
98 #include "gromacs/timing/cyclecounter.h"
99 #include "gromacs/timing/gpu_timing.h"
100 #include "gromacs/timing/wallcycle.h"
101 #include "gromacs/timing/wallcyclereporting.h"
102 #include "gromacs/timing/walltime_accounting.h"
103 #include "gromacs/topology/topology.h"
104 #include "gromacs/utility/arrayref.h"
105 #include "gromacs/utility/basedefinitions.h"
106 #include "gromacs/utility/cstringutil.h"
107 #include "gromacs/utility/exceptions.h"
108 #include "gromacs/utility/fatalerror.h"
109 #include "gromacs/utility/gmxassert.h"
110 #include "gromacs/utility/gmxmpi.h"
111 #include "gromacs/utility/logger.h"
112 #include "gromacs/utility/smalloc.h"
113 #include "gromacs/utility/strconvert.h"
114 #include "gromacs/utility/sysinfo.h"
116 // TODO: this environment variable allows us to verify before release
117 // that on less common architectures the total cost of polling is not larger than
118 // a blocking wait (so polling does not introduce overhead when the static
119 // PME-first ordering would suffice).
120 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
123 static void sum_forces(rvec f[], gmx::ArrayRef<const gmx::RVec> forceToAdd)
125 const int end = forceToAdd.size();
127 int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
128 #pragma omp parallel for num_threads(nt) schedule(static)
129 for (int i = 0; i < end; i++)
131 rvec_inc(f[i], forceToAdd[i]);
135 static void calc_virial(int start, int homenr, const rvec x[], const rvec f[],
136 tensor vir_part, const t_graph *graph, const matrix box,
137 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
139 /* The short-range virial from surrounding boxes */
140 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
141 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
143 /* Calculate partial virial, for local atoms only, based on short range.
144 * Total virial is computed in global_stat, called from do_md
146 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
147 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
151 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
155 static void pull_potential_wrapper(const t_commrec *cr,
156 const t_inputrec *ir,
157 const matrix box, gmx::ArrayRef<const gmx::RVec> x,
158 gmx::ForceWithVirial *force,
159 const t_mdatoms *mdatoms,
160 gmx_enerdata_t *enerd,
163 gmx_wallcycle_t wcycle)
168 /* Calculate the center of mass forces, this requires communication,
169 * which is why pull_potential is called close to other communication.
171 wallcycle_start(wcycle, ewcPULLPOT);
172 set_pbc(&pbc, ir->ePBC, box);
174 enerd->term[F_COM_PULL] +=
175 pull_potential(ir->pull_work, mdatoms, &pbc,
176 cr, t, lambda[efptRESTRAINT], as_rvec_array(x.data()), force, &dvdl);
177 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
178 wallcycle_stop(wcycle, ewcPULLPOT);
181 static void pme_receive_force_ener(const t_commrec *cr,
182 gmx::ForceWithVirial *forceWithVirial,
183 gmx_enerdata_t *enerd,
184 gmx_wallcycle_t wcycle)
186 real e_q, e_lj, dvdl_q, dvdl_lj;
187 float cycles_ppdpme, cycles_seppme;
189 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
190 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
192 /* In case of node-splitting, the PP nodes receive the long-range
193 * forces, virial and energy from the PME nodes here.
195 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
198 gmx_pme_receive_f(cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
200 enerd->term[F_COUL_RECIP] += e_q;
201 enerd->term[F_LJ_RECIP] += e_lj;
202 enerd->dvdl_lin[efptCOUL] += dvdl_q;
203 enerd->dvdl_lin[efptVDW] += dvdl_lj;
207 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
209 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
212 static void print_large_forces(FILE *fp,
220 real force2Tolerance = gmx::square(forceTolerance);
221 gmx::index numNonFinite = 0;
222 for (int i = 0; i < md->homenr; i++)
224 real force2 = norm2(f[i]);
225 bool nonFinite = !std::isfinite(force2);
226 if (force2 >= force2Tolerance || nonFinite)
228 fprintf(fp, "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
230 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
237 if (numNonFinite > 0)
239 /* Note that with MPI this fatal call on one rank might interrupt
240 * the printing on other ranks. But we can only avoid that with
241 * an expensive MPI barrier that we would need at each step.
243 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
247 static void post_process_forces(const t_commrec *cr,
250 gmx_wallcycle_t wcycle,
251 const gmx_localtop_t *top,
255 gmx::ForceWithVirial *forceWithVirial,
257 const t_mdatoms *mdatoms,
258 const t_graph *graph,
259 const t_forcerec *fr,
260 const gmx_vsite_t *vsite,
263 if (fr->haveDirectVirialContributions)
265 rvec *fDirectVir = as_rvec_array(forceWithVirial->force_.data());
269 /* Spread the mesh force on virtual sites to the other particles...
270 * This is parallellized. MPI communication is performed
271 * if the constructing atoms aren't local.
273 matrix virial = { { 0 } };
274 spread_vsite_f(vsite, x, fDirectVir, nullptr,
275 (flags & GMX_FORCE_VIRIAL) != 0, virial,
277 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
278 forceWithVirial->addVirialContribution(virial);
281 if (flags & GMX_FORCE_VIRIAL)
283 /* Now add the forces, this is local */
284 sum_forces(f, forceWithVirial->force_);
286 /* Add the direct virial contributions */
287 GMX_ASSERT(forceWithVirial->computeVirial_, "forceWithVirial should request virial computation when we request the virial");
288 m_add(vir_force, forceWithVirial->getVirial(), vir_force);
292 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
297 if (fr->print_force >= 0)
299 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
303 static void do_nb_verlet(t_forcerec *fr,
304 const interaction_const_t *ic,
305 gmx_enerdata_t *enerd,
307 const Nbnxm::InteractionLocality ilocality,
311 gmx_wallcycle_t wcycle)
313 if (!(flags & GMX_FORCE_NONBONDED))
315 /* skip non-bonded calculation */
319 nonbonded_verlet_t *nbv = fr->nbv.get();
321 /* GPU kernel launch overhead is already timed separately */
322 if (fr->cutoff_scheme != ecutsVERLET)
324 gmx_incons("Invalid cut-off scheme passed!");
329 /* When dynamic pair-list pruning is requested, we need to prune
330 * at nstlistPrune steps.
332 if (nbv->isDynamicPruningStepCpu(step))
334 /* Prune the pair-list beyond fr->ic->rlistPrune using
335 * the current coordinates of the atoms.
337 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
338 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
339 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
342 wallcycle_sub_start(wcycle, ewcsNONBONDED);
345 nbv->dispatchNonbondedKernel(ilocality, *ic, flags, clearF, fr, enerd, nrnb);
349 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
353 static inline void clear_rvecs_omp(int n, rvec v[])
355 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, n);
357 /* Note that we would like to avoid this conditional by putting it
358 * into the omp pragma instead, but then we still take the full
359 * omp parallel for overhead (at least with gcc5).
363 for (int i = 0; i < n; i++)
370 #pragma omp parallel for num_threads(nth) schedule(static)
371 for (int i = 0; i < n; i++)
378 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
380 * \param groupOptions Group options, containing T-coupling options
382 static real averageKineticEnergyEstimate(const t_grpopts &groupOptions)
384 real nrdfCoupled = 0;
385 real nrdfUncoupled = 0;
386 real kineticEnergy = 0;
387 for (int g = 0; g < groupOptions.ngtc; g++)
389 if (groupOptions.tau_t[g] >= 0)
391 nrdfCoupled += groupOptions.nrdf[g];
392 kineticEnergy += groupOptions.nrdf[g]*0.5*groupOptions.ref_t[g]*BOLTZ;
396 nrdfUncoupled += groupOptions.nrdf[g];
400 /* This conditional with > also catches nrdf=0 */
401 if (nrdfCoupled > nrdfUncoupled)
403 return kineticEnergy*(nrdfCoupled + nrdfUncoupled)/nrdfCoupled;
411 /*! \brief This routine checks that the potential energy is finite.
413 * Always checks that the potential energy is finite. If step equals
414 * inputrec.init_step also checks that the magnitude of the potential energy
415 * is reasonable. Terminates with a fatal error when a check fails.
416 * Note that passing this check does not guarantee finite forces,
417 * since those use slightly different arithmetics. But in most cases
418 * there is just a narrow coordinate range where forces are not finite
419 * and energies are finite.
421 * \param[in] step The step number, used for checking and printing
422 * \param[in] enerd The energy data; the non-bonded group energies need to be added to enerd.term[F_EPOT] before calling this routine
423 * \param[in] inputrec The input record
425 static void checkPotentialEnergyValidity(int64_t step,
426 const gmx_enerdata_t &enerd,
427 const t_inputrec &inputrec)
429 /* Threshold valid for comparing absolute potential energy against
430 * the kinetic energy. Normally one should not consider absolute
431 * potential energy values, but with a factor of one million
432 * we should never get false positives.
434 constexpr real c_thresholdFactor = 1e6;
436 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
437 real averageKineticEnergy = 0;
438 /* We only check for large potential energy at the initial step,
439 * because that is by far the most likely step for this too occur
440 * and because computing the average kinetic energy is not free.
441 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
442 * before they become NaN.
444 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
446 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
449 if (energyIsNotFinite || (averageKineticEnergy > 0 &&
450 enerd.term[F_EPOT] > c_thresholdFactor*averageKineticEnergy))
452 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.",
455 energyIsNotFinite ? "not finite" : "extremely high",
457 enerd.term[F_COUL_SR],
458 energyIsNotFinite ? "non-finite" : "very high",
459 energyIsNotFinite ? " or Nan" : "");
463 /*! \brief Return true if there are special forces computed this step.
465 * The conditionals exactly correspond to those in computeSpecialForces().
468 haveSpecialForces(const t_inputrec *inputrec,
469 ForceProviders *forceProviders,
473 const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
476 ((computeForces && forceProviders->hasForceProvider()) || // forceProviders
477 (inputrec->bPull && pull_have_potential(inputrec->pull_work)) || // pull
478 inputrec->bRot || // enforced rotation
479 (ed != nullptr) || // flooding
480 (inputrec->bIMD && computeForces)); // IMD
483 /*! \brief Compute forces and/or energies for special algorithms
485 * The intention is to collect all calls to algorithms that compute
486 * forces on local atoms only and that do not contribute to the local
487 * virial sum (but add their virial contribution separately).
488 * Eventually these should likely all become ForceProviders.
489 * Within this function the intention is to have algorithms that do
490 * global communication at the end, so global barriers within the MD loop
491 * are as close together as possible.
493 * \param[in] fplog The log file
494 * \param[in] cr The communication record
495 * \param[in] inputrec The input record
496 * \param[in] awh The Awh module (nullptr if none in use).
497 * \param[in] enforcedRotation Enforced rotation module.
498 * \param[in] imdSession The IMD session
499 * \param[in] step The current MD step
500 * \param[in] t The current time
501 * \param[in,out] wcycle Wallcycle accounting struct
502 * \param[in,out] forceProviders Pointer to a list of force providers
503 * \param[in] box The unit cell
504 * \param[in] x The coordinates
505 * \param[in] mdatoms Per atom properties
506 * \param[in] lambda Array of free-energy lambda values
507 * \param[in] forceFlags Flags that tell whether we should compute forces/energies/virial
508 * \param[in,out] forceWithVirial Force and virial buffers
509 * \param[in,out] enerd Energy buffer
510 * \param[in,out] ed Essential dynamics pointer
511 * \param[in] bNS Tells if we did neighbor searching this step, used for ED sampling
513 * \todo Remove bNS, which is used incorrectly.
514 * \todo Convert all other algorithms called here to ForceProviders.
517 computeSpecialForces(FILE *fplog,
519 const t_inputrec *inputrec,
521 gmx_enfrot *enforcedRotation,
522 gmx::ImdSession *imdSession,
525 gmx_wallcycle_t wcycle,
526 ForceProviders *forceProviders,
528 gmx::ArrayRef<const gmx::RVec> x,
529 const t_mdatoms *mdatoms,
532 gmx::ForceWithVirial *forceWithVirial,
533 gmx_enerdata_t *enerd,
537 const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
539 /* NOTE: Currently all ForceProviders only provide forces.
540 * When they also provide energies, remove this conditional.
544 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
545 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
547 /* Collect forces from modules */
548 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
551 if (inputrec->bPull && pull_have_potential(inputrec->pull_work))
553 pull_potential_wrapper(cr, inputrec, box, x,
555 mdatoms, enerd, lambda, t,
560 enerd->term[F_COM_PULL] +=
561 awh->applyBiasForcesAndUpdateBias(inputrec->ePBC, *mdatoms, box,
563 t, step, wcycle, fplog);
567 rvec *f = as_rvec_array(forceWithVirial->force_.data());
569 /* Add the forces from enforced rotation potentials (if any) */
572 wallcycle_start(wcycle, ewcROTadd);
573 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
574 wallcycle_stop(wcycle, ewcROTadd);
579 /* Note that since init_edsam() is called after the initialization
580 * of forcerec, edsam doesn't request the noVirSum force buffer.
581 * Thus if no other algorithm (e.g. PME) requires it, the forces
582 * here will contribute to the virial.
584 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, bNS);
587 /* Add forces from interactive molecular dynamics (IMD), if any */
588 if (inputrec->bIMD && computeForces)
590 imdSession->applyForces(f);
594 /*! \brief Launch the prepare_step and spread stages of PME GPU.
596 * \param[in] pmedata The PME structure
597 * \param[in] box The box matrix
598 * \param[in] x Coordinate array
599 * \param[in] flags Force flags
600 * \param[in] pmeFlags PME flags
601 * \param[in] wcycle The wallcycle structure
603 static inline void launchPmeGpuSpread(gmx_pme_t *pmedata,
608 gmx_wallcycle_t wcycle)
610 pme_gpu_prepare_computation(pmedata, (flags & GMX_FORCE_DYNAMICBOX) != 0, box, wcycle, pmeFlags);
611 pme_gpu_launch_spread(pmedata, x, wcycle);
614 /*! \brief Launch the FFT and gather stages of PME GPU
616 * This function only implements setting the output forces (no accumulation).
618 * \param[in] pmedata The PME structure
619 * \param[in] wcycle The wallcycle structure
621 static void launchPmeGpuFftAndGather(gmx_pme_t *pmedata,
622 gmx_wallcycle_t wcycle)
624 pme_gpu_launch_complex_transforms(pmedata, wcycle);
625 pme_gpu_launch_gather(pmedata, wcycle, PmeForceOutputHandling::Set);
629 * Polling wait for either of the PME or nonbonded GPU tasks.
631 * Instead of a static order in waiting for GPU tasks, this function
632 * polls checking which of the two tasks completes first, and does the
633 * associated force buffer reduction overlapped with the other task.
634 * By doing that, unlike static scheduling order, it can always overlap
635 * one of the reductions, regardless of the GPU task completion order.
637 * \param[in] nbv Nonbonded verlet structure
638 * \param[in,out] pmedata PME module data
639 * \param[in,out] force Force array to reduce task outputs into.
640 * \param[in,out] forceWithVirial Force and virial buffers
641 * \param[in,out] fshift Shift force output vector results are reduced into
642 * \param[in,out] enerd Energy data structure results are reduced into
643 * \param[in] flags Force flags
644 * \param[in] pmeFlags PME flags
645 * \param[in] haveOtherWork Tells whether there is other work than non-bonded in the stream(s)
646 * \param[in] wcycle The wallcycle structure
648 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t *nbv,
650 gmx::ArrayRefWithPadding<gmx::RVec> *force,
651 gmx::ForceWithVirial *forceWithVirial,
653 gmx_enerdata_t *enerd,
657 gmx_wallcycle_t wcycle)
659 bool isPmeGpuDone = false;
660 bool isNbGpuDone = false;
663 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
665 while (!isPmeGpuDone || !isNbGpuDone)
669 GpuTaskCompletion completionType = (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
670 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, pmeFlags, wcycle, forceWithVirial, enerd, completionType);
675 GpuTaskCompletion completionType = (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
676 wallcycle_start_nocount(wcycle, ewcWAIT_GPU_NB_L);
677 isNbGpuDone = Nbnxm::gpu_try_finish_task(nbv->gpu_nbv,
679 Nbnxm::AtomLocality::Local,
681 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
682 fshift, completionType);
683 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
684 // To get the call count right, when the task finished we
685 // issue a start/stop.
686 // TODO: move the ewcWAIT_GPU_NB_L cycle counting into nbnxn_gpu_try_finish_task()
687 // and ewcNB_XF_BUF_OPS counting into nbnxn_atomdata_add_nbat_f_to_f().
690 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
691 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
693 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
694 as_rvec_array(force->unpaddedArrayRef().data()), wcycle);
700 /*! \brief Hack structure with force ouput buffers for do_force for the home atoms for this domain */
704 ForceOutputs(rvec *f, gmx::ForceWithVirial const forceWithVirial) :
706 forceWithVirial(forceWithVirial) {}
708 //! Force output buffer used by legacy modules (without SIMD padding)
710 //! Force with direct virial contribution (if there are any; without SIMD padding)
711 gmx::ForceWithVirial forceWithVirial;
714 /*! \brief Set up the different force buffers; also does clearing.
716 * \param[in] fr force record pointer
717 * \param[in] inputrec input record
718 * \param[in] force force array
719 * \param[in] bDoForces True if force are computed this step
720 * \param[in] doVirial True if virial is computed this step
721 * \param[out] wcycle wallcycle recording structure
723 * \returns Cleared force output structure
726 setupForceOutputs(const t_forcerec *fr,
727 const t_inputrec &inputrec,
728 gmx::ArrayRefWithPadding<gmx::RVec> force,
729 const bool bDoForces,
731 gmx_wallcycle_t wcycle)
733 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
735 /* Temporary solution until all routines take PaddedRVecVector */
736 rvec *const f = as_rvec_array(force.unpaddedArrayRef().data());
739 /* Clear the short- and long-range forces */
740 clear_rvecs_omp(fr->natoms_force_constr, f);
743 /* If we need to compute the virial, we might need a separate
744 * force buffer for algorithms for which the virial is calculated
745 * directly, such as PME. Otherwise, forceWithVirial uses the
746 * the same force (f in legacy calls) buffer as other algorithms.
748 const bool useSeparateForceWithVirialBuffer = (bDoForces && (doVirial && fr->haveDirectVirialContributions));
751 /* forceWithVirial uses the local atom range only */
752 gmx::ForceWithVirial forceWithVirial (useSeparateForceWithVirialBuffer ?
753 *fr->forceBufferForDirectVirialContributions : force.unpaddedArrayRef(),
756 if (useSeparateForceWithVirialBuffer)
758 /* TODO: update comment
759 * We only compute forces on local atoms. Note that vsites can
760 * spread to non-local atoms, but that part of the buffer is
761 * cleared separately in the vsite spreading code.
763 clear_rvecs_omp(forceWithVirial.force_.size(), as_rvec_array(forceWithVirial.force_.data()));
766 if (inputrec.bPull && pull_have_constraint(inputrec.pull_work))
768 clear_pull_forces(inputrec.pull_work);
771 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
773 return ForceOutputs(f, forceWithVirial);
777 /*! \brief Set up flags that indicate what type of work is there to compute.
779 * Currently we only update it at search steps,
780 * but some properties may change more frequently (e.g. virial/non-virial step),
781 * so when including those either the frequency of update (per-step) or the scope
782 * of a flag will change (i.e. a set of flags for nstlist steps).
786 setupForceWorkload(gmx::PpForceWorkload *forceWork,
787 const t_inputrec *inputrec,
788 const t_forcerec *fr,
795 forceWork->haveSpecialForces = haveSpecialForces(inputrec, fr->forceProviders, forceFlags, ed);
796 forceWork->haveCpuBondedWork = haveCpuBondeds(*fr);
797 forceWork->haveGpuBondedWork = ((fr->gpuBonded != nullptr) && fr->gpuBonded->haveInteractions());
798 forceWork->haveRestraintsWork = havePositionRestraints(idef, *fcd);
799 forceWork->haveCpuListedForceWork = haveCpuListedForces(*fr, idef, *fcd);
802 void do_force(FILE *fplog,
804 const gmx_multisim_t *ms,
805 const t_inputrec *inputrec,
807 gmx_enfrot *enforcedRotation,
808 gmx::ImdSession *imdSession,
811 gmx_wallcycle_t wcycle,
812 const gmx_localtop_t *top,
814 gmx::ArrayRefWithPadding<gmx::RVec> x, //NOLINT(performance-unnecessary-value-param)
816 gmx::ArrayRefWithPadding<gmx::RVec> force, //NOLINT(performance-unnecessary-value-param)
818 const t_mdatoms *mdatoms,
819 gmx_enerdata_t *enerd,
821 gmx::ArrayRef<real> lambda,
824 gmx::PpForceWorkload *ppForceWorkload,
825 const gmx_vsite_t *vsite,
830 const DDBalanceRegionHandler &ddBalanceRegionHandler)
834 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
835 gmx_bool bDoForces, bUseGPU, bUseOrEmulGPU;
836 nonbonded_verlet_t *nbv = fr->nbv.get();
837 interaction_const_t *ic = fr->ic;
839 /* modify force flag if not doing nonbonded */
842 flags &= ~GMX_FORCE_NONBONDED;
844 bStateChanged = ((flags & GMX_FORCE_STATECHANGED) != 0);
845 bNS = ((flags & GMX_FORCE_NS) != 0);
846 bFillGrid = (bNS && bStateChanged);
847 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
848 bDoForces = ((flags & GMX_FORCE_FORCES) != 0);
849 bUseGPU = fr->nbv->useGpu();
850 bUseOrEmulGPU = bUseGPU || fr->nbv->emulateGpu();
852 const auto pmeRunMode = fr->pmedata ? pme_run_mode(fr->pmedata) : PmeRunMode::CPU;
853 // TODO slim this conditional down - inputrec and duty checks should mean the same in proper code!
854 const bool useGpuPme = EEL_PME(fr->ic->eeltype) && thisRankHasDuty(cr, DUTY_PME) &&
855 ((pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed));
856 const int pmeFlags = GMX_PME_SPREAD | GMX_PME_SOLVE |
857 ((flags & GMX_FORCE_VIRIAL) ? GMX_PME_CALC_ENER_VIR : 0) |
858 ((flags & GMX_FORCE_ENERGY) ? GMX_PME_CALC_ENER_VIR : 0) |
859 ((flags & GMX_FORCE_FORCES) ? GMX_PME_CALC_F : 0);
861 /* At a search step we need to start the first balancing region
862 * somewhere early inside the step after communication during domain
863 * decomposition (and not during the previous step as usual).
867 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
871 const int homenr = mdatoms->homenr;
873 clear_mat(vir_force);
877 update_forcerec(fr, box);
879 if (inputrecNeedMutot(inputrec))
881 /* Calculate total (local) dipole moment in a temporary common array.
882 * This makes it possible to sum them over nodes faster.
884 calc_mu(start, homenr,
885 x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
890 if (fr->ePBC != epbcNONE)
892 /* Compute shift vectors every step,
893 * because of pressure coupling or box deformation!
895 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
897 calc_shifts(box, fr->shift_vec);
902 put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr), gmx_omp_nthreads_get(emntDefault));
903 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
905 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
907 unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
911 nbnxn_atomdata_copy_shiftvec((flags & GMX_FORCE_DYNAMICBOX) != 0,
912 fr->shift_vec, nbv->nbat.get());
915 if (!thisRankHasDuty(cr, DUTY_PME))
917 /* Send particle coordinates to the pme nodes.
918 * Since this is only implemented for domain decomposition
919 * and domain decomposition does not use the graph,
920 * we do not need to worry about shifting.
922 gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
923 lambda[efptCOUL], lambda[efptVDW],
924 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)) != 0,
931 launchPmeGpuSpread(fr->pmedata, box, as_rvec_array(x.unpaddedArrayRef().data()), flags, pmeFlags, wcycle);
934 /* do gridding for pair search */
937 if (graph && bStateChanged)
939 /* Calculate intramolecular shift vectors to make molecules whole */
940 mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
944 // - vzero is constant, do we need to pass it?
945 // - box_diag should be passed directly to nbnxn_put_on_grid
951 box_diag[XX] = box[XX][XX];
952 box_diag[YY] = box[YY][YY];
953 box_diag[ZZ] = box[ZZ][ZZ];
955 wallcycle_start(wcycle, ewcNS);
956 if (!DOMAINDECOMP(cr))
958 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
959 nbnxn_put_on_grid(nbv, box,
961 nullptr, 0, mdatoms->homenr, -1,
962 fr->cginfo, x.unpaddedArrayRef(),
964 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
968 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
969 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd),
970 fr->cginfo, x.unpaddedArrayRef());
971 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
974 nbv->setAtomProperties(*mdatoms, *fr->cginfo);
976 wallcycle_stop(wcycle, ewcNS);
978 /* initialize the GPU nbnxm atom data and bonded data structures */
981 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
983 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
984 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
985 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
989 /* Now we put all atoms on the grid, we can assign bonded
990 * interactions to the GPU, where the grid order is
991 * needed. Also the xq, f and fshift device buffers have
992 * been reallocated if needed, so the bonded code can
993 * learn about them. */
994 // TODO the xq, f, and fshift buffers are now shared
995 // resources, so they should be maintained by a
996 // higher-level object than the nb module.
997 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
999 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1000 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1001 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1003 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1006 // Need to run after the GPU-offload bonded interaction lists
1007 // are set up to be able to determine whether there is bonded work.
1008 setupForceWorkload(ppForceWorkload,
1017 /* do local pair search */
1020 // TODO: fuse this branch with the above bNS block
1021 wallcycle_start_nocount(wcycle, ewcNS);
1022 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1023 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1024 nbv->constructPairlist(Nbnxm::InteractionLocality::Local,
1025 &top->excls, step, nrnb);
1026 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1027 wallcycle_stop(wcycle, ewcNS);
1031 nbv->setCoordinates(Nbnxm::AtomLocality::Local, false,
1032 x.unpaddedArrayRef(), wcycle);
1037 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1039 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1041 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1042 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1043 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1044 Nbnxm::AtomLocality::Local,
1045 ppForceWorkload->haveGpuBondedWork);
1046 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1048 // bonded work not split into separate local and non-local, so with DD
1049 // we can only launch the kernel after non-local coordinates have been received.
1050 if (ppForceWorkload->haveGpuBondedWork && !havePPDomainDecomposition(cr))
1052 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1053 fr->gpuBonded->launchKernels(fr, flags, box);
1054 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1057 /* launch local nonbonded work on GPU */
1058 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1059 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFNo,
1060 step, nrnb, wcycle);
1061 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1062 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1067 // In PME GPU and mixed mode we launch FFT / gather after the
1068 // X copy/transform to allow overlap as well as after the GPU NB
1069 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1070 // the nonbonded kernel.
1071 launchPmeGpuFftAndGather(fr->pmedata, wcycle);
1074 /* Communicate coordinates and sum dipole if necessary +
1075 do non-local pair search */
1076 if (havePPDomainDecomposition(cr))
1080 // TODO: fuse this branch with the above large bNS block
1081 wallcycle_start_nocount(wcycle, ewcNS);
1082 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1083 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1084 nbv->constructPairlist(Nbnxm::InteractionLocality::NonLocal,
1085 &top->excls, step, nrnb);
1086 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1087 wallcycle_stop(wcycle, ewcNS);
1091 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1093 nbv->setCoordinates(Nbnxm::AtomLocality::NonLocal, false,
1094 x.unpaddedArrayRef(), wcycle);
1099 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1101 /* launch non-local nonbonded tasks on GPU */
1102 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1103 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1104 Nbnxm::AtomLocality::NonLocal,
1105 ppForceWorkload->haveGpuBondedWork);
1106 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1108 if (ppForceWorkload->haveGpuBondedWork)
1110 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1111 fr->gpuBonded->launchKernels(fr, flags, box);
1112 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1115 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1116 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1117 step, nrnb, wcycle);
1118 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1120 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1126 /* launch D2H copy-back F */
1127 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1128 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1129 if (havePPDomainDecomposition(cr))
1131 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1132 flags, Nbnxm::AtomLocality::NonLocal, ppForceWorkload->haveGpuBondedWork);
1134 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1135 flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork);
1136 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1138 if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1140 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1141 fr->gpuBonded->launchEnergyTransfer();
1142 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1144 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1147 if (bStateChanged && inputrecNeedMutot(inputrec))
1151 gmx_sumd(2*DIM, mu, cr);
1153 ddBalanceRegionHandler.reopenRegionCpu();
1156 for (i = 0; i < 2; i++)
1158 for (j = 0; j < DIM; j++)
1160 fr->mu_tot[i][j] = mu[i*DIM + j];
1164 if (fr->efep == efepNO)
1166 copy_rvec(fr->mu_tot[0], mu_tot);
1170 for (j = 0; j < DIM; j++)
1173 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1174 lambda[efptCOUL]*fr->mu_tot[1][j];
1178 /* Reset energies */
1179 reset_enerdata(enerd);
1180 clear_rvecs(SHIFTS, fr->fshift);
1182 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1184 wallcycle_start(wcycle, ewcPPDURINGPME);
1185 dd_force_flop_start(cr->dd, nrnb);
1190 wallcycle_start(wcycle, ewcROT);
1191 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, bNS);
1192 wallcycle_stop(wcycle, ewcROT);
1195 /* Start the force cycle counter.
1196 * Note that a different counter is used for dynamic load balancing.
1198 wallcycle_start(wcycle, ewcFORCE);
1200 // set up and clear force outputs
1201 struct ForceOutputs forceOut = setupForceOutputs(fr, *inputrec, force, bDoForces, ((flags & GMX_FORCE_VIRIAL) != 0), wcycle);
1203 /* We calculate the non-bonded forces, when done on the CPU, here.
1204 * We do this before calling do_force_lowlevel, because in that
1205 * function, the listed forces are calculated before PME, which
1206 * does communication. With this order, non-bonded and listed
1207 * force calculation imbalance can be balanced out by the domain
1208 * decomposition load balancing.
1213 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFYes,
1214 step, nrnb, wcycle);
1217 if (fr->efep != efepNO)
1219 /* Calculate the local and non-local free energy interactions here.
1220 * Happens here on the CPU both with and without GPU.
1222 wallcycle_sub_start(wcycle, ewcsNONBONDED);
1223 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::Local,
1224 fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1225 inputrec->fepvals, lambda.data(),
1226 enerd, flags, nrnb);
1228 if (havePPDomainDecomposition(cr))
1230 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::NonLocal,
1231 fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1232 inputrec->fepvals, lambda.data(),
1233 enerd, flags, nrnb);
1235 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
1240 if (havePPDomainDecomposition(cr))
1242 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1243 step, nrnb, wcycle);
1246 /* Add all the non-bonded force to the normal force array.
1247 * This can be split into a local and a non-local part when overlapping
1248 * communication with calculation with domain decomposition.
1250 wallcycle_stop(wcycle, ewcFORCE);
1252 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::All, forceOut.f, wcycle);
1254 wallcycle_start_nocount(wcycle, ewcFORCE);
1256 /* If there are multiple fshift output buffers we need to reduce them */
1257 if (flags & GMX_FORCE_VIRIAL)
1259 /* This is not in a subcounter because it takes a
1260 negligible and constant-sized amount of time */
1261 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->nbat.get(),
1266 /* update QMMMrec, if necessary */
1269 update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1272 /* Compute the bonded and non-bonded energies and optionally forces */
1273 do_force_lowlevel(fr, inputrec, &(top->idef),
1274 cr, ms, nrnb, wcycle, mdatoms,
1275 as_rvec_array(x.unpaddedArrayRef().data()), hist, forceOut.f, &forceOut.forceWithVirial, enerd, fcd,
1276 box, inputrec->fepvals, lambda.data(), graph, &(top->excls), fr->mu_tot,
1278 ddBalanceRegionHandler);
1280 wallcycle_stop(wcycle, ewcFORCE);
1282 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1283 imdSession, step, t, wcycle,
1284 fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda.data(),
1285 flags, &forceOut.forceWithVirial, enerd,
1288 // Will store the amount of cycles spent waiting for the GPU that
1289 // will be later used in the DLB accounting.
1290 float cycles_wait_gpu = 0;
1293 /* wait for non-local forces (or calculate in emulation mode) */
1294 if (havePPDomainDecomposition(cr))
1298 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1299 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1300 flags, Nbnxm::AtomLocality::NonLocal,
1301 ppForceWorkload->haveGpuBondedWork,
1302 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1304 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1308 wallcycle_start_nocount(wcycle, ewcFORCE);
1309 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFYes,
1310 step, nrnb, wcycle);
1311 wallcycle_stop(wcycle, ewcFORCE);
1314 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::NonLocal,
1315 forceOut.f, wcycle);
1319 if (havePPDomainDecomposition(cr))
1321 /* We are done with the CPU compute.
1322 * We will now communicate the non-local forces.
1323 * If we use a GPU this will overlap with GPU work, so in that case
1324 * we do not close the DD force balancing region here.
1326 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1330 dd_move_f(cr->dd, force.unpaddedArrayRef(), fr->fshift, wcycle);
1334 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1335 // an alternating wait/reduction scheme.
1336 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPme && bUseGPU && !DOMAINDECOMP(cr));
1337 if (alternateGpuWait)
1339 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &force, &forceOut.forceWithVirial, fr->fshift, enerd,
1340 flags, pmeFlags, ppForceWorkload->haveGpuBondedWork, wcycle);
1343 if (!alternateGpuWait && useGpuPme)
1345 pme_gpu_wait_and_reduce(fr->pmedata, pmeFlags, wcycle, &forceOut.forceWithVirial, enerd);
1348 /* Wait for local GPU NB outputs on the non-alternating wait path */
1349 if (!alternateGpuWait && bUseGPU)
1351 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1352 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1353 * but even with a step of 0.1 ms the difference is less than 1%
1356 const float gpuWaitApiOverheadMargin = 2e6f; /* cycles */
1358 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1359 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1360 flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork,
1361 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1363 float cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1365 if (ddBalanceRegionHandler.useBalancingRegion())
1367 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1368 if (bDoForces && cycles_tmp <= gpuWaitApiOverheadMargin)
1370 /* We measured few cycles, it could be that the kernel
1371 * and transfer finished earlier and there was no actual
1372 * wait time, only API call overhead.
1373 * Then the actual time could be anywhere between 0 and
1374 * cycles_wait_est. We will use half of cycles_wait_est.
1376 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1378 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1382 if (fr->nbv->emulateGpu())
1384 // NOTE: emulation kernel is not included in the balancing region,
1385 // but emulation mode does not target performance anyway
1386 wallcycle_start_nocount(wcycle, ewcFORCE);
1387 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local,
1388 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1389 step, nrnb, wcycle);
1390 wallcycle_stop(wcycle, ewcFORCE);
1395 pme_gpu_reinit_computation(fr->pmedata, wcycle);
1400 /* now clear the GPU outputs while we finish the step on the CPU */
1401 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1402 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1403 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, flags);
1405 if (nbv->isDynamicPruningStepGpu(step))
1407 nbv->dispatchPruneKernelGpu(step);
1409 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1410 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1413 if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1415 wallcycle_start(wcycle, ewcWAIT_GPU_BONDED);
1416 // in principle this should be included in the DD balancing region,
1417 // but generally it is infrequent so we'll omit it for the sake of
1419 fr->gpuBonded->accumulateEnergyTerms(enerd);
1420 wallcycle_stop(wcycle, ewcWAIT_GPU_BONDED);
1422 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1423 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1424 fr->gpuBonded->clearEnergies();
1425 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1426 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1429 /* Do the nonbonded GPU (or emulation) force buffer reduction
1430 * on the non-alternating path. */
1431 if (bUseOrEmulGPU && !alternateGpuWait)
1433 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
1434 forceOut.f, wcycle);
1436 if (DOMAINDECOMP(cr))
1438 dd_force_flop_stop(cr->dd, nrnb);
1443 /* If we have NoVirSum forces, but we do not calculate the virial,
1444 * we sum fr->f_novirsum=forceOut.f later.
1446 if (vsite && !(fr->haveDirectVirialContributions && !(flags & GMX_FORCE_VIRIAL)))
1448 spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, fr->fshift, FALSE, nullptr, nrnb,
1449 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1452 if (flags & GMX_FORCE_VIRIAL)
1454 /* Calculation of the virial must be done after vsites! */
1455 calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f,
1456 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1460 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1462 /* In case of node-splitting, the PP nodes receive the long-range
1463 * forces, virial and energy from the PME nodes here.
1465 pme_receive_force_ener(cr, &forceOut.forceWithVirial, enerd, wcycle);
1470 post_process_forces(cr, step, nrnb, wcycle,
1471 top, box, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, &forceOut.forceWithVirial,
1472 vir_force, mdatoms, graph, fr, vsite,
1476 if (flags & GMX_FORCE_ENERGY)
1478 /* Sum the potential energy terms from group contributions */
1479 sum_epot(&(enerd->grpp), enerd->term);
1481 if (!EI_TPI(inputrec->eI))
1483 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1487 /* In case we don't have constraints and are using GPUs, the next balancing
1488 * region starts here.
1489 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1490 * virial calculation and COM pulling, is not thus not included in
1491 * the balance timing, which is ok as most tasks do communication.
1493 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);