2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013-2019,2020, 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/gpuhaloexchange.h"
53 #include "gromacs/domdec/partition.h"
54 #include "gromacs/essentialdynamics/edsam.h"
55 #include "gromacs/ewald/pme.h"
56 #include "gromacs/ewald/pme_pp.h"
57 #include "gromacs/ewald/pme_pp_comm_gpu.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
60 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
61 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
62 #include "gromacs/gpu_utils/gpu_utils.h"
63 #include "gromacs/imd/imd.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/force_flags.h"
80 #include "gromacs/mdlib/forcerec.h"
81 #include "gromacs/mdlib/gmx_omp_nthreads.h"
82 #include "gromacs/mdlib/qmmm.h"
83 #include "gromacs/mdlib/update.h"
84 #include "gromacs/mdlib/vsite.h"
85 #include "gromacs/mdtypes/commrec.h"
86 #include "gromacs/mdtypes/enerdata.h"
87 #include "gromacs/mdtypes/forceoutput.h"
88 #include "gromacs/mdtypes/forcerec.h"
89 #include "gromacs/mdtypes/iforceprovider.h"
90 #include "gromacs/mdtypes/inputrec.h"
91 #include "gromacs/mdtypes/md_enums.h"
92 #include "gromacs/mdtypes/mdatom.h"
93 #include "gromacs/mdtypes/simulation_workload.h"
94 #include "gromacs/mdtypes/state.h"
95 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
96 #include "gromacs/nbnxm/gpu_data_mgmt.h"
97 #include "gromacs/nbnxm/nbnxm.h"
98 #include "gromacs/nbnxm/nbnxm_gpu.h"
99 #include "gromacs/pbcutil/ishift.h"
100 #include "gromacs/pbcutil/mshift.h"
101 #include "gromacs/pbcutil/pbc.h"
102 #include "gromacs/pulling/pull.h"
103 #include "gromacs/pulling/pull_rotation.h"
104 #include "gromacs/timing/cyclecounter.h"
105 #include "gromacs/timing/gpu_timing.h"
106 #include "gromacs/timing/wallcycle.h"
107 #include "gromacs/timing/wallcyclereporting.h"
108 #include "gromacs/timing/walltime_accounting.h"
109 #include "gromacs/topology/topology.h"
110 #include "gromacs/utility/arrayref.h"
111 #include "gromacs/utility/basedefinitions.h"
112 #include "gromacs/utility/cstringutil.h"
113 #include "gromacs/utility/exceptions.h"
114 #include "gromacs/utility/fatalerror.h"
115 #include "gromacs/utility/fixedcapacityvector.h"
116 #include "gromacs/utility/gmxassert.h"
117 #include "gromacs/utility/gmxmpi.h"
118 #include "gromacs/utility/logger.h"
119 #include "gromacs/utility/smalloc.h"
120 #include "gromacs/utility/strconvert.h"
121 #include "gromacs/utility/sysinfo.h"
123 using gmx::AtomLocality;
124 using gmx::DomainLifetimeWorkload;
125 using gmx::ForceOutputs;
126 using gmx::InteractionLocality;
127 using gmx::SimulationWorkload;
128 using gmx::StepWorkload;
130 // TODO: this environment variable allows us to verify before release
131 // that on less common architectures the total cost of polling is not larger than
132 // a blocking wait (so polling does not introduce overhead when the static
133 // PME-first ordering would suffice).
134 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
136 static void sum_forces(rvec f[], gmx::ArrayRef<const gmx::RVec> forceToAdd)
138 const int end = forceToAdd.size();
140 int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
141 #pragma omp parallel for num_threads(nt) schedule(static)
142 for (int i = 0; i < end; i++)
144 rvec_inc(f[i], forceToAdd[i]);
148 static void calc_virial(int start,
151 const gmx::ForceWithShiftForces& forceWithShiftForces,
153 const t_graph* graph,
156 const t_forcerec* fr,
159 /* The short-range virial from surrounding boxes */
160 const rvec* fshift = as_rvec_array(forceWithShiftForces.shiftForces().data());
161 calc_vir(SHIFTS, fr->shift_vec, fshift, vir_part, pbcType == PbcType::Screw, box);
162 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
164 /* Calculate partial virial, for local atoms only, based on short range.
165 * Total virial is computed in global_stat, called from do_md
167 const rvec* f = as_rvec_array(forceWithShiftForces.force().data());
168 f_calc_vir(start, start + homenr, x, f, vir_part, graph, box);
169 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
173 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
177 static void pull_potential_wrapper(const t_commrec* cr,
178 const t_inputrec* ir,
180 gmx::ArrayRef<const gmx::RVec> x,
181 gmx::ForceWithVirial* force,
182 const t_mdatoms* mdatoms,
183 gmx_enerdata_t* enerd,
187 gmx_wallcycle_t wcycle)
192 /* Calculate the center of mass forces, this requires communication,
193 * which is why pull_potential is called close to other communication.
195 wallcycle_start(wcycle, ewcPULLPOT);
196 set_pbc(&pbc, ir->pbcType, box);
198 enerd->term[F_COM_PULL] += pull_potential(pull_work, mdatoms, &pbc, cr, t, lambda[efptRESTRAINT],
199 as_rvec_array(x.data()), force, &dvdl);
200 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
201 for (auto& dhdl : enerd->dhdlLambda)
205 wallcycle_stop(wcycle, ewcPULLPOT);
208 static void pme_receive_force_ener(t_forcerec* fr,
210 gmx::ForceWithVirial* forceWithVirial,
211 gmx_enerdata_t* enerd,
212 bool useGpuPmePpComms,
213 bool receivePmeForceToGpu,
214 gmx_wallcycle_t wcycle)
216 real e_q, e_lj, dvdl_q, dvdl_lj;
217 float cycles_ppdpme, cycles_seppme;
219 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
220 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
222 /* In case of node-splitting, the PP nodes receive the long-range
223 * forces, virial and energy from the PME nodes here.
225 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
228 gmx_pme_receive_f(fr->pmePpCommGpu.get(), cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
229 useGpuPmePpComms, receivePmeForceToGpu, &cycles_seppme);
230 enerd->term[F_COUL_RECIP] += e_q;
231 enerd->term[F_LJ_RECIP] += e_lj;
232 enerd->dvdl_lin[efptCOUL] += dvdl_q;
233 enerd->dvdl_lin[efptVDW] += dvdl_lj;
235 for (auto& dhdl : enerd->dhdlLambda)
237 dhdl += dvdl_q + dvdl_lj;
242 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
244 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
247 static void print_large_forces(FILE* fp,
255 real force2Tolerance = gmx::square(forceTolerance);
256 gmx::index numNonFinite = 0;
257 for (int i = 0; i < md->homenr; i++)
259 real force2 = norm2(f[i]);
260 bool nonFinite = !std::isfinite(force2);
261 if (force2 >= force2Tolerance || nonFinite)
263 fprintf(fp, "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n", step,
264 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
271 if (numNonFinite > 0)
273 /* Note that with MPI this fatal call on one rank might interrupt
274 * the printing on other ranks. But we can only avoid that with
275 * an expensive MPI barrier that we would need at each step.
277 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
281 static void post_process_forces(const t_commrec* cr,
284 gmx_wallcycle_t wcycle,
285 const gmx_localtop_t* top,
288 ForceOutputs* forceOutputs,
290 const t_mdatoms* mdatoms,
291 const t_graph* graph,
292 const t_forcerec* fr,
293 const gmx_vsite_t* vsite,
294 const StepWorkload& stepWork)
296 rvec* f = as_rvec_array(forceOutputs->forceWithShiftForces().force().data());
298 if (fr->haveDirectVirialContributions)
300 auto& forceWithVirial = forceOutputs->forceWithVirial();
301 rvec* fDirectVir = as_rvec_array(forceWithVirial.force_.data());
305 /* Spread the mesh force on virtual sites to the other particles...
306 * This is parallellized. MPI communication is performed
307 * if the constructing atoms aren't local.
309 matrix virial = { { 0 } };
310 spread_vsite_f(vsite, x, fDirectVir, nullptr, stepWork.computeVirial, virial, nrnb,
311 top->idef, fr->pbcType, fr->bMolPBC, graph, box, cr, wcycle);
312 forceWithVirial.addVirialContribution(virial);
315 if (stepWork.computeVirial)
317 /* Now add the forces, this is local */
318 sum_forces(f, forceWithVirial.force_);
320 /* Add the direct virial contributions */
322 forceWithVirial.computeVirial_,
323 "forceWithVirial should request virial computation when we request the virial");
324 m_add(vir_force, forceWithVirial.getVirial(), vir_force);
328 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
333 if (fr->print_force >= 0)
335 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
339 static void do_nb_verlet(t_forcerec* fr,
340 const interaction_const_t* ic,
341 gmx_enerdata_t* enerd,
342 const StepWorkload& stepWork,
343 const InteractionLocality ilocality,
347 gmx_wallcycle_t wcycle)
349 if (!stepWork.computeNonbondedForces)
351 /* skip non-bonded calculation */
355 nonbonded_verlet_t* nbv = fr->nbv.get();
357 /* GPU kernel launch overhead is already timed separately */
358 if (fr->cutoff_scheme != ecutsVERLET)
360 gmx_incons("Invalid cut-off scheme passed!");
365 /* When dynamic pair-list pruning is requested, we need to prune
366 * at nstlistPrune steps.
368 if (nbv->isDynamicPruningStepCpu(step))
370 /* Prune the pair-list beyond fr->ic->rlistPrune using
371 * the current coordinates of the atoms.
373 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
374 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
375 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
379 nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
382 static inline void clear_rvecs_omp(int n, rvec v[])
384 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, n);
386 /* Note that we would like to avoid this conditional by putting it
387 * into the omp pragma instead, but then we still take the full
388 * omp parallel for overhead (at least with gcc5).
392 for (int i = 0; i < n; i++)
399 #pragma omp parallel for num_threads(nth) schedule(static)
400 for (int i = 0; i < n; i++)
407 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
409 * \param groupOptions Group options, containing T-coupling options
411 static real averageKineticEnergyEstimate(const t_grpopts& groupOptions)
413 real nrdfCoupled = 0;
414 real nrdfUncoupled = 0;
415 real kineticEnergy = 0;
416 for (int g = 0; g < groupOptions.ngtc; g++)
418 if (groupOptions.tau_t[g] >= 0)
420 nrdfCoupled += groupOptions.nrdf[g];
421 kineticEnergy += groupOptions.nrdf[g] * 0.5 * groupOptions.ref_t[g] * BOLTZ;
425 nrdfUncoupled += groupOptions.nrdf[g];
429 /* This conditional with > also catches nrdf=0 */
430 if (nrdfCoupled > nrdfUncoupled)
432 return kineticEnergy * (nrdfCoupled + nrdfUncoupled) / nrdfCoupled;
440 /*! \brief This routine checks that the potential energy is finite.
442 * Always checks that the potential energy is finite. If step equals
443 * inputrec.init_step also checks that the magnitude of the potential energy
444 * is reasonable. Terminates with a fatal error when a check fails.
445 * Note that passing this check does not guarantee finite forces,
446 * since those use slightly different arithmetics. But in most cases
447 * there is just a narrow coordinate range where forces are not finite
448 * and energies are finite.
450 * \param[in] step The step number, used for checking and printing
451 * \param[in] enerd The energy data; the non-bonded group energies need to be added to
452 * enerd.term[F_EPOT] before calling this routine \param[in] inputrec The input record
454 static void checkPotentialEnergyValidity(int64_t step, const gmx_enerdata_t& enerd, const t_inputrec& inputrec)
456 /* Threshold valid for comparing absolute potential energy against
457 * the kinetic energy. Normally one should not consider absolute
458 * potential energy values, but with a factor of one million
459 * we should never get false positives.
461 constexpr real c_thresholdFactor = 1e6;
463 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
464 real averageKineticEnergy = 0;
465 /* We only check for large potential energy at the initial step,
466 * because that is by far the most likely step for this too occur
467 * and because computing the average kinetic energy is not free.
468 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
469 * before they become NaN.
471 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
473 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
476 if (energyIsNotFinite
477 || (averageKineticEnergy > 0 && enerd.term[F_EPOT] > c_thresholdFactor * averageKineticEnergy))
482 ": The total potential energy is %g, which is %s. The LJ and electrostatic "
483 "contributions to the energy are %g and %g, respectively. A %s potential energy "
484 "can be caused by overlapping interactions in bonded interactions or very large%s "
485 "coordinate values. Usually this is caused by a badly- or non-equilibrated initial "
486 "configuration, incorrect interactions or parameters in the topology.",
487 step, enerd.term[F_EPOT], energyIsNotFinite ? "not finite" : "extremely high",
488 enerd.term[F_LJ], enerd.term[F_COUL_SR],
489 energyIsNotFinite ? "non-finite" : "very high", energyIsNotFinite ? " or Nan" : "");
493 /*! \brief Return true if there are special forces computed this step.
495 * The conditionals exactly correspond to those in computeSpecialForces().
497 static bool haveSpecialForces(const t_inputrec& inputrec,
498 const gmx::ForceProviders& forceProviders,
499 const pull_t* pull_work,
500 const bool computeForces,
504 return ((computeForces && forceProviders.hasForceProvider()) || // forceProviders
505 (inputrec.bPull && pull_have_potential(pull_work)) || // pull
506 inputrec.bRot || // enforced rotation
507 (ed != nullptr) || // flooding
508 (inputrec.bIMD && computeForces)); // IMD
511 /*! \brief Compute forces and/or energies for special algorithms
513 * The intention is to collect all calls to algorithms that compute
514 * forces on local atoms only and that do not contribute to the local
515 * virial sum (but add their virial contribution separately).
516 * Eventually these should likely all become ForceProviders.
517 * Within this function the intention is to have algorithms that do
518 * global communication at the end, so global barriers within the MD loop
519 * are as close together as possible.
521 * \param[in] fplog The log file
522 * \param[in] cr The communication record
523 * \param[in] inputrec The input record
524 * \param[in] awh The Awh module (nullptr if none in use).
525 * \param[in] enforcedRotation Enforced rotation module.
526 * \param[in] imdSession The IMD session
527 * \param[in] pull_work The pull work structure.
528 * \param[in] step The current MD step
529 * \param[in] t The current time
530 * \param[in,out] wcycle Wallcycle accounting struct
531 * \param[in,out] forceProviders Pointer to a list of force providers
532 * \param[in] box The unit cell
533 * \param[in] x The coordinates
534 * \param[in] mdatoms Per atom properties
535 * \param[in] lambda Array of free-energy lambda values
536 * \param[in] stepWork Step schedule flags
537 * \param[in,out] forceWithVirial Force and virial buffers
538 * \param[in,out] enerd Energy buffer
539 * \param[in,out] ed Essential dynamics pointer
540 * \param[in] didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
542 * \todo Remove didNeighborSearch, which is used incorrectly.
543 * \todo Convert all other algorithms called here to ForceProviders.
545 static void computeSpecialForces(FILE* fplog,
547 const t_inputrec* inputrec,
549 gmx_enfrot* enforcedRotation,
550 gmx::ImdSession* imdSession,
554 gmx_wallcycle_t wcycle,
555 gmx::ForceProviders* forceProviders,
557 gmx::ArrayRef<const gmx::RVec> x,
558 const t_mdatoms* mdatoms,
560 const StepWorkload& stepWork,
561 gmx::ForceWithVirial* forceWithVirial,
562 gmx_enerdata_t* enerd,
564 bool didNeighborSearch)
566 /* NOTE: Currently all ForceProviders only provide forces.
567 * When they also provide energies, remove this conditional.
569 if (stepWork.computeForces)
571 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
572 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
574 /* Collect forces from modules */
575 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
578 if (inputrec->bPull && pull_have_potential(pull_work))
580 pull_potential_wrapper(cr, inputrec, box, x, forceWithVirial, mdatoms, enerd, pull_work,
585 enerd->term[F_COM_PULL] += awh->applyBiasForcesAndUpdateBias(
586 inputrec->pbcType, *mdatoms, box, forceWithVirial, t, step, wcycle, fplog);
590 rvec* f = as_rvec_array(forceWithVirial->force_.data());
592 /* Add the forces from enforced rotation potentials (if any) */
595 wallcycle_start(wcycle, ewcROTadd);
596 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
597 wallcycle_stop(wcycle, ewcROTadd);
602 /* Note that since init_edsam() is called after the initialization
603 * of forcerec, edsam doesn't request the noVirSum force buffer.
604 * Thus if no other algorithm (e.g. PME) requires it, the forces
605 * here will contribute to the virial.
607 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
610 /* Add forces from interactive molecular dynamics (IMD), if any */
611 if (inputrec->bIMD && stepWork.computeForces)
613 imdSession->applyForces(f);
617 /*! \brief Makes PME flags from StepWorkload data.
619 * \param[in] stepWork Step schedule flags
622 static int makePmeFlags(const StepWorkload& stepWork)
624 return GMX_PME_SPREAD | GMX_PME_SOLVE | (stepWork.computeVirial ? GMX_PME_CALC_ENER_VIR : 0)
625 | (stepWork.computeEnergy ? GMX_PME_CALC_ENER_VIR : 0)
626 | (stepWork.computeForces ? GMX_PME_CALC_F : 0);
629 /*! \brief Launch the prepare_step and spread stages of PME GPU.
631 * \param[in] pmedata The PME structure
632 * \param[in] box The box matrix
633 * \param[in] stepWork Step schedule flags
634 * \param[in] pmeFlags PME flags
635 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in
636 * the device memory. \param[in] wcycle The wallcycle structure
638 static inline void launchPmeGpuSpread(gmx_pme_t* pmedata,
640 const StepWorkload& stepWork,
642 GpuEventSynchronizer* xReadyOnDevice,
643 gmx_wallcycle_t wcycle)
645 pme_gpu_prepare_computation(pmedata, stepWork.haveDynamicBox, box, wcycle, pmeFlags,
646 stepWork.useGpuPmeFReduction);
647 pme_gpu_launch_spread(pmedata, xReadyOnDevice, wcycle);
650 /*! \brief Launch the FFT and gather stages of PME GPU
652 * This function only implements setting the output forces (no accumulation).
654 * \param[in] pmedata The PME structure
655 * \param[in] wcycle The wallcycle structure
657 static void launchPmeGpuFftAndGather(gmx_pme_t* pmedata, gmx_wallcycle_t wcycle)
659 pme_gpu_launch_complex_transforms(pmedata, wcycle);
660 pme_gpu_launch_gather(pmedata, wcycle);
664 * Polling wait for either of the PME or nonbonded GPU tasks.
666 * Instead of a static order in waiting for GPU tasks, this function
667 * polls checking which of the two tasks completes first, and does the
668 * associated force buffer reduction overlapped with the other task.
669 * By doing that, unlike static scheduling order, it can always overlap
670 * one of the reductions, regardless of the GPU task completion order.
672 * \param[in] nbv Nonbonded verlet structure
673 * \param[in,out] pmedata PME module data
674 * \param[in,out] forceOutputs Output buffer for the forces and virial
675 * \param[in,out] enerd Energy data structure results are reduced into
676 * \param[in] stepWork Step schedule flags
677 * \param[in] pmeFlags PME flags
678 * \param[in] wcycle The wallcycle structure
680 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t* nbv,
682 gmx::ForceOutputs* forceOutputs,
683 gmx_enerdata_t* enerd,
684 const StepWorkload& stepWork,
686 gmx_wallcycle_t wcycle)
688 bool isPmeGpuDone = false;
689 bool isNbGpuDone = false;
692 gmx::ForceWithShiftForces& forceWithShiftForces = forceOutputs->forceWithShiftForces();
693 gmx::ForceWithVirial& forceWithVirial = forceOutputs->forceWithVirial();
695 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
697 while (!isPmeGpuDone || !isNbGpuDone)
701 GpuTaskCompletion completionType =
702 (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
703 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, pmeFlags, wcycle, &forceWithVirial,
704 enerd, completionType);
709 GpuTaskCompletion completionType =
710 (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
711 isNbGpuDone = Nbnxm::gpu_try_finish_task(
712 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
713 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(),
714 completionType, wcycle);
718 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShiftForces.force());
724 /*! \brief Set up the different force buffers; also does clearing.
726 * \param[in] fr force record pointer
727 * \param[in] pull_work The pull work object.
728 * \param[in] inputrec input record
729 * \param[in] force force array
730 * \param[in] stepWork Step schedule flags
731 * \param[out] wcycle wallcycle recording structure
733 * \returns Cleared force output structure
735 static ForceOutputs setupForceOutputs(t_forcerec* fr,
737 const t_inputrec& inputrec,
738 gmx::ArrayRefWithPadding<gmx::RVec> force,
739 const StepWorkload& stepWork,
740 gmx_wallcycle_t wcycle)
742 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
744 /* NOTE: We assume fr->shiftForces is all zeros here */
745 gmx::ForceWithShiftForces forceWithShiftForces(force, stepWork.computeVirial, fr->shiftForces);
747 if (stepWork.computeForces)
749 /* Clear the short- and long-range forces */
750 clear_rvecs_omp(fr->natoms_force_constr, as_rvec_array(forceWithShiftForces.force().data()));
753 /* If we need to compute the virial, we might need a separate
754 * force buffer for algorithms for which the virial is calculated
755 * directly, such as PME. Otherwise, forceWithVirial uses the
756 * the same force (f in legacy calls) buffer as other algorithms.
758 const bool useSeparateForceWithVirialBuffer =
759 (stepWork.computeForces && (stepWork.computeVirial && fr->haveDirectVirialContributions));
760 /* forceWithVirial uses the local atom range only */
761 gmx::ForceWithVirial forceWithVirial(useSeparateForceWithVirialBuffer ? fr->forceBufferForDirectVirialContributions
762 : force.unpaddedArrayRef(),
763 stepWork.computeVirial);
765 if (useSeparateForceWithVirialBuffer)
767 /* TODO: update comment
768 * We only compute forces on local atoms. Note that vsites can
769 * spread to non-local atoms, but that part of the buffer is
770 * cleared separately in the vsite spreading code.
772 clear_rvecs_omp(forceWithVirial.force_.size(), as_rvec_array(forceWithVirial.force_.data()));
775 if (inputrec.bPull && pull_have_constraint(pull_work))
777 clear_pull_forces(pull_work);
780 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
782 return ForceOutputs(forceWithShiftForces, forceWithVirial);
786 /*! \brief Set up flags that have the lifetime of the domain indicating what type of work is there to compute.
788 static DomainLifetimeWorkload setupDomainLifetimeWorkload(const t_inputrec& inputrec,
789 const t_forcerec& fr,
790 const pull_t* pull_work,
792 const InteractionDefinitions& idef,
794 const t_mdatoms& mdatoms,
795 const SimulationWorkload& simulationWork,
796 const StepWorkload& stepWork)
798 DomainLifetimeWorkload domainWork;
799 // Note that haveSpecialForces is constant over the whole run
800 domainWork.haveSpecialForces =
801 haveSpecialForces(inputrec, *fr.forceProviders, pull_work, stepWork.computeForces, ed);
802 domainWork.haveCpuBondedWork = haveCpuBondeds(fr);
803 domainWork.haveGpuBondedWork = ((fr.gpuBonded != nullptr) && fr.gpuBonded->haveInteractions());
804 domainWork.haveRestraintsWork = haveRestraints(idef, fcd);
805 domainWork.haveCpuListedForceWork = haveCpuListedForces(fr, idef, fcd);
806 // Note that haveFreeEnergyWork is constant over the whole run
807 domainWork.haveFreeEnergyWork = (fr.efep != efepNO && mdatoms.nPerturbed != 0);
808 // We assume we have local force work if there are CPU
809 // force tasks including PME or nonbondeds.
810 domainWork.haveCpuLocalForceWork =
811 domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork
812 || domainWork.haveFreeEnergyWork || simulationWork.useCpuNonbonded || simulationWork.useCpuPme
813 || simulationWork.haveEwaldSurfaceContribution || inputrec.nwall > 0;
818 /*! \brief Set up force flag stuct from the force bitmask.
820 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
821 * \param[in] isNonbondedOn Global override, if false forces to turn off all nonbonded calculation.
822 * \param[in] simulationWork Simulation workload description.
823 * \param[in] rankHasPmeDuty If this rank computes PME.
825 * \returns New Stepworkload description.
827 static StepWorkload setupStepWorkload(const int legacyFlags,
828 const bool isNonbondedOn,
829 const SimulationWorkload& simulationWork,
830 const bool rankHasPmeDuty)
833 flags.stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
834 flags.haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
835 flags.doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
836 flags.computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
837 flags.computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
838 flags.computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
839 flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
840 flags.computeNonbondedForces = ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && isNonbondedOn;
841 flags.computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
843 if (simulationWork.useGpuBufferOps)
845 GMX_ASSERT(simulationWork.useGpuNonbonded,
846 "Can only offload buffer ops if nonbonded computation is also offloaded");
848 flags.useGpuXBufferOps = simulationWork.useGpuBufferOps;
849 // on virial steps the CPU reduction path is taken
850 flags.useGpuFBufferOps = simulationWork.useGpuBufferOps && !flags.computeVirial;
851 flags.useGpuPmeFReduction = flags.useGpuFBufferOps
852 && (simulationWork.useGpuPme
853 && (rankHasPmeDuty || simulationWork.useGpuPmePpCommunication));
859 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
861 * TODO: eliminate \p useGpuPmeOnThisRank when this is
862 * incorporated in DomainLifetimeWorkload.
864 static void launchGpuEndOfStepTasks(nonbonded_verlet_t* nbv,
865 gmx::GpuBonded* gpuBonded,
867 gmx_enerdata_t* enerd,
868 const gmx::MdrunScheduleWorkload& runScheduleWork,
869 bool useGpuPmeOnThisRank,
871 gmx_wallcycle_t wcycle)
873 if (runScheduleWork.simulationWork.useGpuNonbonded)
875 /* Launch pruning before buffer clearing because the API overhead of the
876 * clear kernel launches can leave the GPU idle while it could be running
879 if (nbv->isDynamicPruningStepGpu(step))
881 nbv->dispatchPruneKernelGpu(step);
884 /* now clear the GPU outputs while we finish the step on the CPU */
885 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
886 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
887 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, runScheduleWork.stepWork.computeVirial);
888 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
889 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
892 if (useGpuPmeOnThisRank)
894 pme_gpu_reinit_computation(pmedata, wcycle);
897 if (runScheduleWork.domainWork.haveGpuBondedWork && runScheduleWork.stepWork.computeEnergy)
899 // in principle this should be included in the DD balancing region,
900 // but generally it is infrequent so we'll omit it for the sake of
902 gpuBonded->waitAccumulateEnergyTerms(enerd);
904 gpuBonded->clearEnergies();
909 void do_force(FILE* fplog,
911 const gmx_multisim_t* ms,
912 const t_inputrec* inputrec,
914 gmx_enfrot* enforcedRotation,
915 gmx::ImdSession* imdSession,
919 gmx_wallcycle_t wcycle,
920 const gmx_localtop_t* top,
922 gmx::ArrayRefWithPadding<gmx::RVec> x,
924 gmx::ArrayRefWithPadding<gmx::RVec> force,
926 const t_mdatoms* mdatoms,
927 gmx_enerdata_t* enerd,
929 gmx::ArrayRef<real> lambda,
932 gmx::MdrunScheduleWorkload* runScheduleWork,
933 const gmx_vsite_t* vsite,
938 const DDBalanceRegionHandler& ddBalanceRegionHandler)
942 nonbonded_verlet_t* nbv = fr->nbv.get();
943 interaction_const_t* ic = fr->ic;
944 gmx::StatePropagatorDataGpu* stateGpu = fr->stateGpu;
946 const SimulationWorkload& simulationWork = runScheduleWork->simulationWork;
949 runScheduleWork->stepWork = setupStepWorkload(legacyFlags, fr->bNonbonded, simulationWork,
950 thisRankHasDuty(cr, DUTY_PME));
951 const StepWorkload& stepWork = runScheduleWork->stepWork;
954 const bool useGpuPmeOnThisRank = simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME);
955 const int pmeFlags = makePmeFlags(stepWork);
957 /* At a search step we need to start the first balancing region
958 * somewhere early inside the step after communication during domain
959 * decomposition (and not during the previous step as usual).
961 if (stepWork.doNeighborSearch)
963 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
967 const int homenr = mdatoms->homenr;
969 clear_mat(vir_force);
971 if (stepWork.stateChanged && simulationWork.computeMuTot)
973 /* Calculate total (local) dipole moment in a temporary common array.
974 * This makes it possible to sum them over nodes faster.
976 calc_mu(start, homenr, x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB,
977 mdatoms->nChargePerturbed, mu, mu + DIM);
980 if (fr->pbcType != PbcType::No)
982 /* Compute shift vectors every step,
983 * because of pressure coupling or box deformation!
985 if (stepWork.haveDynamicBox && stepWork.stateChanged)
987 calc_shifts(box, fr->shift_vec);
990 const bool fillGrid = (stepWork.doNeighborSearch && stepWork.stateChanged);
991 const bool calcCGCM = (fillGrid && !DOMAINDECOMP(cr));
994 put_atoms_in_box_omp(fr->pbcType, box, x.unpaddedArrayRef().subArray(0, homenr),
995 gmx_omp_nthreads_get(emntDefault));
996 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
998 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1000 unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
1004 nbnxn_atomdata_copy_shiftvec(stepWork.haveDynamicBox, fr->shift_vec, nbv->nbat.get());
1007 const bool pmeSendCoordinatesFromGpu =
1008 simulationWork.useGpuPmePpCommunication && !(stepWork.doNeighborSearch);
1009 const bool reinitGpuPmePpComms =
1010 simulationWork.useGpuPmePpCommunication && (stepWork.doNeighborSearch);
1013 const auto localXReadyOnDevice = (stateGpu != nullptr)
1014 ? stateGpu->getCoordinatesReadyOnDeviceEvent(
1015 AtomLocality::Local, simulationWork, stepWork)
1019 // If coordinates are to be sent to PME task from CPU memory, perform that send here.
1020 // Otherwise the send will occur after H2D coordinate transfer.
1021 if (!thisRankHasDuty(cr, DUTY_PME) && !pmeSendCoordinatesFromGpu)
1023 /* Send particle coordinates to the pme nodes.
1024 * Since this is only implemented for domain decomposition
1025 * and domain decomposition does not use the graph,
1026 * we do not need to worry about shifting.
1028 if (!stepWork.doNeighborSearch && simulationWork.useGpuUpdate)
1030 GMX_RELEASE_ASSERT(false,
1031 "GPU update and separate PME ranks are only supported with GPU "
1032 "direct communication!");
1033 // TODO: when this code-path becomes supported add:
1034 // stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1037 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1038 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1039 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1040 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1042 #endif /* GMX_MPI */
1044 // Coordinates on the device are needed if PME or BufferOps are offloaded.
1045 // The local coordinates can be copied right away.
1046 // NOTE: Consider moving this copy to right after they are updated and constrained,
1047 // if the later is not offloaded.
1048 if (useGpuPmeOnThisRank || stepWork.useGpuXBufferOps)
1050 if (stepWork.doNeighborSearch)
1052 // TODO refactor this to do_md, after partitioning.
1053 stateGpu->reinit(mdatoms->homenr,
1054 cr->dd != nullptr ? dd_numAtomsZones(*cr->dd) : mdatoms->homenr);
1055 if (useGpuPmeOnThisRank)
1057 // TODO: This should be moved into PME setup function ( pme_gpu_prepare_computation(...) )
1058 pme_gpu_set_device_x(fr->pmedata, stateGpu->getCoordinates());
1061 // We need to copy coordinates when:
1062 // 1. Update is not offloaded
1063 // 2. The buffers were reinitialized on search step
1064 if (!simulationWork.useGpuUpdate || stepWork.doNeighborSearch)
1066 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1067 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1071 // TODO Update this comment when introducing SimulationWorkload
1073 // The conditions for gpuHaloExchange e.g. using GPU buffer
1074 // operations were checked before construction, so here we can
1075 // just use it and assert upon any conditions.
1076 const bool ddUsesGpuDirectCommunication =
1077 ((cr->dd != nullptr) && (!cr->dd->gpuHaloExchange.empty()));
1078 GMX_ASSERT(!ddUsesGpuDirectCommunication || stepWork.useGpuXBufferOps,
1079 "Must use coordinate buffer ops with GPU halo exchange");
1080 const bool useGpuForcesHaloExchange = ddUsesGpuDirectCommunication && stepWork.useGpuFBufferOps;
1082 // Copy coordinate from the GPU if update is on the GPU and there
1083 // are forces to be computed on the CPU, or for the computation of
1084 // virial, or if host-side data will be transferred from this task
1085 // to a remote task for halo exchange or PME-PP communication. At
1086 // search steps the current coordinates are already on the host,
1087 // hence copy is not needed.
1088 const bool haveHostPmePpComms =
1089 !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication;
1090 const bool haveHostHaloExchangeComms = havePPDomainDecomposition(cr) && !ddUsesGpuDirectCommunication;
1092 bool gmx_used_in_debug haveCopiedXFromGpu = false;
1093 if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
1094 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
1095 || haveHostPmePpComms || haveHostHaloExchangeComms))
1097 GMX_ASSERT(stateGpu != nullptr, "stateGpu should not be null");
1098 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
1099 haveCopiedXFromGpu = true;
1103 // If coordinates are to be sent to PME task from GPU memory, perform that send here.
1104 // Otherwise the send will occur before the H2D coordinate transfer.
1105 if (pmeSendCoordinatesFromGpu)
1107 /* Send particle coordinates to the pme nodes.
1108 * Since this is only implemented for domain decomposition
1109 * and domain decomposition does not use the graph,
1110 * we do not need to worry about shifting.
1112 gmx_pme_send_coordinates(fr, cr, box, as_rvec_array(x.unpaddedArrayRef().data()), lambda[efptCOUL],
1113 lambda[efptVDW], (stepWork.computeVirial || stepWork.computeEnergy),
1114 step, simulationWork.useGpuPmePpCommunication, reinitGpuPmePpComms,
1115 pmeSendCoordinatesFromGpu, localXReadyOnDevice, wcycle);
1117 #endif /* GMX_MPI */
1119 if (useGpuPmeOnThisRank)
1121 launchPmeGpuSpread(fr->pmedata, box, stepWork, pmeFlags, localXReadyOnDevice, wcycle);
1124 /* do gridding for pair search */
1125 if (stepWork.doNeighborSearch)
1127 if (graph && stepWork.stateChanged)
1129 /* Calculate intramolecular shift vectors to make molecules whole */
1130 mk_mshift(fplog, graph, fr->pbcType, box, as_rvec_array(x.unpaddedArrayRef().data()));
1134 // - vzero is constant, do we need to pass it?
1135 // - box_diag should be passed directly to nbnxn_put_on_grid
1141 box_diag[XX] = box[XX][XX];
1142 box_diag[YY] = box[YY][YY];
1143 box_diag[ZZ] = box[ZZ][ZZ];
1145 wallcycle_start(wcycle, ewcNS);
1146 if (!DOMAINDECOMP(cr))
1148 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1149 nbnxn_put_on_grid(nbv, box, 0, vzero, box_diag, nullptr, { 0, mdatoms->homenr }, -1,
1150 fr->cginfo, x.unpaddedArrayRef(), 0, nullptr);
1151 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1155 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1156 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd), fr->cginfo, x.unpaddedArrayRef());
1157 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1160 nbv->setAtomProperties(*mdatoms, fr->cginfo);
1162 wallcycle_stop(wcycle, ewcNS);
1164 /* initialize the GPU nbnxm atom data and bonded data structures */
1165 if (simulationWork.useGpuNonbonded)
1167 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1169 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1170 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1171 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1175 /* Now we put all atoms on the grid, we can assign bonded
1176 * interactions to the GPU, where the grid order is
1177 * needed. Also the xq, f and fshift device buffers have
1178 * been reallocated if needed, so the bonded code can
1179 * learn about them. */
1180 // TODO the xq, f, and fshift buffers are now shared
1181 // resources, so they should be maintained by a
1182 // higher-level object than the nb module.
1183 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(
1184 nbv->getGridIndices(), top->idef, Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1185 Nbnxm::gpu_get_f(nbv->gpu_nbv), Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1187 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1190 // Need to run after the GPU-offload bonded interaction lists
1191 // are set up to be able to determine whether there is bonded work.
1192 runScheduleWork->domainWork = setupDomainLifetimeWorkload(
1193 *inputrec, *fr, pull_work, ed, top->idef, *fcd, *mdatoms, simulationWork, stepWork);
1195 wallcycle_start_nocount(wcycle, ewcNS);
1196 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1197 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1198 nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
1200 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::Local);
1202 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1203 wallcycle_stop(wcycle, ewcNS);
1205 if (stepWork.useGpuXBufferOps)
1207 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1209 // For force buffer ops, we use the below conditon rather than
1210 // useGpuFBufferOps to ensure that init is performed even if this
1211 // NS step is also a virial step (on which f buf ops are deactivated).
1212 if (simulationWork.useGpuBufferOps && simulationWork.useGpuNonbonded && (GMX_GPU == GMX_GPU_CUDA))
1214 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1215 nbv->atomdata_init_add_nbat_f_to_f_gpu(stateGpu->fReducedOnDevice());
1218 else if (!EI_TPI(inputrec->eI))
1220 if (stepWork.useGpuXBufferOps)
1222 GMX_ASSERT(stateGpu, "stateGpu should be valid when buffer ops are offloaded");
1223 nbv->convertCoordinatesGpu(AtomLocality::Local, false, stateGpu->getCoordinates(),
1224 localXReadyOnDevice);
1228 if (simulationWork.useGpuUpdate)
1230 GMX_ASSERT(stateGpu, "need a valid stateGpu object");
1231 GMX_ASSERT(haveCopiedXFromGpu,
1232 "a wait should only be triggered if copy has been scheduled");
1233 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1235 nbv->convertCoordinates(AtomLocality::Local, false, x.unpaddedArrayRef());
1239 const gmx::DomainLifetimeWorkload& domainWork = runScheduleWork->domainWork;
1241 if (simulationWork.useGpuNonbonded)
1243 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1245 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1247 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1248 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1249 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1251 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::Local);
1253 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1254 // with X buffer ops offloaded to the GPU on all but the search steps
1256 // bonded work not split into separate local and non-local, so with DD
1257 // we can only launch the kernel after non-local coordinates have been received.
1258 if (domainWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1260 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1261 fr->gpuBonded->launchKernel(fr, stepWork, box);
1262 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1265 /* launch local nonbonded work on GPU */
1266 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1267 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFNo, step, nrnb, wcycle);
1268 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1269 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1272 if (useGpuPmeOnThisRank)
1274 // In PME GPU and mixed mode we launch FFT / gather after the
1275 // X copy/transform to allow overlap as well as after the GPU NB
1276 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1277 // the nonbonded kernel.
1278 launchPmeGpuFftAndGather(fr->pmedata, wcycle);
1281 /* Communicate coordinates and sum dipole if necessary +
1282 do non-local pair search */
1283 if (havePPDomainDecomposition(cr))
1285 if (stepWork.doNeighborSearch)
1287 // TODO: fuse this branch with the above large stepWork.doNeighborSearch block
1288 wallcycle_start_nocount(wcycle, ewcNS);
1289 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1290 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1291 nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
1293 nbv->setupGpuShortRangeWork(fr->gpuBonded, InteractionLocality::NonLocal);
1294 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1295 wallcycle_stop(wcycle, ewcNS);
1296 // TODO refactor this GPU halo exchange re-initialisation
1297 // to location in do_md where GPU halo exchange is
1298 // constructed at partitioning, after above stateGpu
1299 // re-initialization has similarly been refactored
1300 if (ddUsesGpuDirectCommunication)
1302 reinitGpuHaloExchange(*cr, stateGpu->getCoordinates(), stateGpu->getForces());
1307 if (ddUsesGpuDirectCommunication)
1309 // The following must be called after local setCoordinates (which records an event
1310 // when the coordinate data has been copied to the device).
1311 communicateGpuHaloCoordinates(*cr, box, localXReadyOnDevice);
1313 if (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork)
1315 // non-local part of coordinate buffer must be copied back to host for CPU work
1316 stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1321 // Note: GPU update + DD without direct communication is not supported,
1322 // a waitCoordinatesReadyOnHost() should be issued if it will be.
1323 GMX_ASSERT(!simulationWork.useGpuUpdate,
1324 "GPU update is not supported with CPU halo exchange");
1325 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1328 if (stepWork.useGpuXBufferOps)
1330 if (!useGpuPmeOnThisRank && !ddUsesGpuDirectCommunication)
1332 stateGpu->copyCoordinatesToGpu(x.unpaddedArrayRef(), AtomLocality::NonLocal);
1334 nbv->convertCoordinatesGpu(AtomLocality::NonLocal, false, stateGpu->getCoordinates(),
1335 stateGpu->getCoordinatesReadyOnDeviceEvent(
1336 AtomLocality::NonLocal, simulationWork, stepWork));
1340 nbv->convertCoordinates(AtomLocality::NonLocal, false, x.unpaddedArrayRef());
1344 if (simulationWork.useGpuNonbonded)
1346 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1348 if (stepWork.doNeighborSearch || !stepWork.useGpuXBufferOps)
1350 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1351 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(), AtomLocality::NonLocal);
1352 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1355 if (domainWork.haveGpuBondedWork)
1357 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1358 fr->gpuBonded->launchKernel(fr, stepWork, box);
1359 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1362 /* launch non-local nonbonded tasks on GPU */
1363 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1364 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1366 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1368 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1372 if (simulationWork.useGpuNonbonded)
1374 /* launch D2H copy-back F */
1375 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1376 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1378 if (havePPDomainDecomposition(cr))
1380 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::NonLocal);
1382 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(), stepWork, AtomLocality::Local);
1383 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1385 if (domainWork.haveGpuBondedWork && stepWork.computeEnergy)
1387 fr->gpuBonded->launchEnergyTransfer();
1389 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1392 if (stepWork.stateChanged && simulationWork.computeMuTot)
1396 gmx_sumd(2 * DIM, mu, cr);
1398 ddBalanceRegionHandler.reopenRegionCpu();
1401 for (i = 0; i < 2; i++)
1403 for (j = 0; j < DIM; j++)
1405 fr->mu_tot[i][j] = mu[i * DIM + j];
1409 if (mdatoms->nChargePerturbed == 0)
1411 copy_rvec(fr->mu_tot[0], mu_tot);
1415 for (j = 0; j < DIM; j++)
1417 mu_tot[j] = (1.0 - lambda[efptCOUL]) * fr->mu_tot[0][j] + lambda[efptCOUL] * fr->mu_tot[1][j];
1421 /* Reset energies */
1422 reset_enerdata(enerd);
1423 /* Clear the shift forces */
1424 // TODO: This should be linked to the shift force buffer in use, or cleared before use instead
1425 for (gmx::RVec& elem : fr->shiftForces)
1427 elem = { 0.0_real, 0.0_real, 0.0_real };
1430 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1432 wallcycle_start(wcycle, ewcPPDURINGPME);
1433 dd_force_flop_start(cr->dd, nrnb);
1436 // For the rest of the CPU tasks that depend on GPU-update produced coordinates,
1437 // this wait ensures that the D2H transfer is complete.
1438 if ((simulationWork.useGpuUpdate)
1439 && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial))
1441 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1446 wallcycle_start(wcycle, ewcROT);
1447 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step,
1448 stepWork.doNeighborSearch);
1449 wallcycle_stop(wcycle, ewcROT);
1452 /* Start the force cycle counter.
1453 * Note that a different counter is used for dynamic load balancing.
1455 wallcycle_start(wcycle, ewcFORCE);
1457 // Set up and clear force outputs.
1458 // We use std::move to keep the compiler happy, it has no effect.
1459 ForceOutputs forceOut =
1460 setupForceOutputs(fr, pull_work, *inputrec, std::move(force), stepWork, wcycle);
1462 /* We calculate the non-bonded forces, when done on the CPU, here.
1463 * We do this before calling do_force_lowlevel, because in that
1464 * function, the listed forces are calculated before PME, which
1465 * does communication. With this order, non-bonded and listed
1466 * force calculation imbalance can be balanced out by the domain
1467 * decomposition load balancing.
1470 const bool useOrEmulateGpuNb = simulationWork.useGpuNonbonded || fr->nbv->emulateGpu();
1472 if (!useOrEmulateGpuNb)
1474 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local, enbvClearFYes, step, nrnb, wcycle);
1477 if (fr->efep != efepNO)
1479 /* Calculate the local and non-local free energy interactions here.
1480 * Happens here on the CPU both with and without GPU.
1482 nbv->dispatchFreeEnergyKernel(InteractionLocality::Local, fr,
1483 as_rvec_array(x.unpaddedArrayRef().data()),
1484 &forceOut.forceWithShiftForces(), *mdatoms, inputrec->fepvals,
1485 lambda.data(), enerd, stepWork, nrnb);
1487 if (havePPDomainDecomposition(cr))
1489 nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal, fr,
1490 as_rvec_array(x.unpaddedArrayRef().data()),
1491 &forceOut.forceWithShiftForces(), *mdatoms,
1492 inputrec->fepvals, lambda.data(), enerd, stepWork, nrnb);
1496 if (!useOrEmulateGpuNb)
1498 if (havePPDomainDecomposition(cr))
1500 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFNo, step,
1504 if (stepWork.computeForces)
1506 /* Add all the non-bonded force to the normal force array.
1507 * This can be split into a local and a non-local part when overlapping
1508 * communication with calculation with domain decomposition.
1510 wallcycle_stop(wcycle, ewcFORCE);
1511 nbv->atomdata_add_nbat_f_to_f(AtomLocality::All, forceOut.forceWithShiftForces().force());
1512 wallcycle_start_nocount(wcycle, ewcFORCE);
1515 /* If there are multiple fshift output buffers we need to reduce them */
1516 if (stepWork.computeVirial)
1518 /* This is not in a subcounter because it takes a
1519 negligible and constant-sized amount of time */
1520 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat,
1521 forceOut.forceWithShiftForces().shiftForces());
1525 /* update QMMMrec, if necessary */
1528 update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1531 // TODO Force flags should include haveFreeEnergyWork for this domain
1532 if (ddUsesGpuDirectCommunication && (domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork))
1534 /* Wait for non-local coordinate data to be copied from device */
1535 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::NonLocal);
1537 /* Compute the bonded and non-bonded energies and optionally forces */
1538 do_force_lowlevel(fr, inputrec, top->idef, cr, ms, nrnb, wcycle, mdatoms, x, hist, &forceOut, enerd,
1539 fcd, box, lambda.data(), graph, fr->mu_tot, stepWork, ddBalanceRegionHandler);
1541 wallcycle_stop(wcycle, ewcFORCE);
1543 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation, imdSession, pull_work, step, t,
1544 wcycle, fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda.data(),
1545 stepWork, &forceOut.forceWithVirial(), enerd, ed, stepWork.doNeighborSearch);
1548 // Will store the amount of cycles spent waiting for the GPU that
1549 // will be later used in the DLB accounting.
1550 float cycles_wait_gpu = 0;
1551 if (useOrEmulateGpuNb)
1553 auto& forceWithShiftForces = forceOut.forceWithShiftForces();
1555 /* wait for non-local forces (or calculate in emulation mode) */
1556 if (havePPDomainDecomposition(cr))
1558 if (simulationWork.useGpuNonbonded)
1560 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(
1561 nbv->gpu_nbv, stepWork, AtomLocality::NonLocal, enerd->grpp.ener[egLJSR].data(),
1562 enerd->grpp.ener[egCOULSR].data(), forceWithShiftForces.shiftForces(), wcycle);
1566 wallcycle_start_nocount(wcycle, ewcFORCE);
1567 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::NonLocal, enbvClearFYes,
1568 step, nrnb, wcycle);
1569 wallcycle_stop(wcycle, ewcFORCE);
1572 if (stepWork.useGpuFBufferOps)
1574 gmx::FixedCapacityVector<GpuEventSynchronizer*, 1> dependencyList;
1576 // TODO: move this into DomainLifetimeWorkload, including the second part of the
1577 // condition The bonded and free energy CPU tasks can have non-local force
1578 // contributions which are a dependency for the GPU force reduction.
1579 bool haveNonLocalForceContribInCpuBuffer =
1580 domainWork.haveCpuBondedWork || domainWork.haveFreeEnergyWork;
1582 if (haveNonLocalForceContribInCpuBuffer)
1584 stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(),
1585 AtomLocality::NonLocal);
1586 dependencyList.push_back(stateGpu->getForcesReadyOnDeviceEvent(
1587 AtomLocality::NonLocal, stepWork.useGpuFBufferOps));
1590 nbv->atomdata_add_nbat_f_to_f_gpu(AtomLocality::NonLocal, stateGpu->getForces(),
1591 pme_gpu_get_device_f(fr->pmedata), dependencyList,
1592 false, haveNonLocalForceContribInCpuBuffer);
1593 if (!useGpuForcesHaloExchange)
1595 // copy from GPU input for dd_move_f()
1596 stateGpu->copyForcesFromGpu(forceOut.forceWithShiftForces().force(),
1597 AtomLocality::NonLocal);
1602 nbv->atomdata_add_nbat_f_to_f(AtomLocality::NonLocal, forceWithShiftForces.force());
1606 if (fr->nbv->emulateGpu() && stepWork.computeVirial)
1608 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat, forceWithShiftForces.shiftForces());
1613 if (havePPDomainDecomposition(cr))
1615 /* We are done with the CPU compute.
1616 * We will now communicate the non-local forces.
1617 * If we use a GPU this will overlap with GPU work, so in that case
1618 * we do not close the DD force balancing region here.
1620 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1622 if (stepWork.computeForces)
1625 if (useGpuForcesHaloExchange)
1627 if (domainWork.haveCpuLocalForceWork)
1629 stateGpu->copyForcesToGpu(forceOut.forceWithShiftForces().force(), AtomLocality::Local);
1631 communicateGpuHaloForces(*cr, domainWork.haveCpuLocalForceWork);
1635 if (stepWork.useGpuFBufferOps)
1637 stateGpu->waitForcesReadyOnHost(AtomLocality::NonLocal);
1639 dd_move_f(cr->dd, &forceOut.forceWithShiftForces(), wcycle);
1644 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1645 // an alternating wait/reduction scheme.
1646 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPmeOnThisRank && simulationWork.useGpuNonbonded
1647 && !DOMAINDECOMP(cr) && !stepWork.useGpuFBufferOps);
1648 if (alternateGpuWait)
1650 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &forceOut, enerd, stepWork, pmeFlags, wcycle);
1653 if (!alternateGpuWait && useGpuPmeOnThisRank)
1655 pme_gpu_wait_and_reduce(fr->pmedata, pmeFlags, wcycle, &forceOut.forceWithVirial(), enerd);
1658 /* Wait for local GPU NB outputs on the non-alternating wait path */
1659 if (!alternateGpuWait && simulationWork.useGpuNonbonded)
1661 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1662 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1663 * but even with a step of 0.1 ms the difference is less than 1%
1666 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
1667 const float waitCycles = Nbnxm::gpu_wait_finish_task(
1668 nbv->gpu_nbv, stepWork, AtomLocality::Local, enerd->grpp.ener[egLJSR].data(),
1669 enerd->grpp.ener[egCOULSR].data(), forceOut.forceWithShiftForces().shiftForces(), wcycle);
1671 if (ddBalanceRegionHandler.useBalancingRegion())
1673 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1674 if (stepWork.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
1676 /* We measured few cycles, it could be that the kernel
1677 * and transfer finished earlier and there was no actual
1678 * wait time, only API call overhead.
1679 * Then the actual time could be anywhere between 0 and
1680 * cycles_wait_est. We will use half of cycles_wait_est.
1682 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1684 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1688 if (fr->nbv->emulateGpu())
1690 // NOTE: emulation kernel is not included in the balancing region,
1691 // but emulation mode does not target performance anyway
1692 wallcycle_start_nocount(wcycle, ewcFORCE);
1693 do_nb_verlet(fr, ic, enerd, stepWork, InteractionLocality::Local,
1694 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes, step, nrnb, wcycle);
1695 wallcycle_stop(wcycle, ewcFORCE);
1698 // If on GPU PME-PP comms or GPU update path, receive forces from PME before GPU buffer ops
1699 // TODO refactor this and unify with below default-path call to the same function
1700 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME)
1701 && (simulationWork.useGpuPmePpCommunication || simulationWork.useGpuUpdate))
1703 /* In case of node-splitting, the PP nodes receive the long-range
1704 * forces, virial and energy from the PME nodes here.
1706 pme_receive_force_ener(fr, cr, &forceOut.forceWithVirial(), enerd,
1707 simulationWork.useGpuPmePpCommunication,
1708 stepWork.useGpuPmeFReduction, wcycle);
1712 /* Do the nonbonded GPU (or emulation) force buffer reduction
1713 * on the non-alternating path. */
1714 if (useOrEmulateGpuNb && !alternateGpuWait)
1716 // TODO simplify the below conditionals. Pass buffer and sync pointers at init stage rather than here. Unify getter fns for sameGPU/otherGPU cases.
1718 stepWork.useGpuPmeFReduction
1719 ? (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_device_f(fr->pmedata)
1720 : // PME force buffer on same GPU
1721 fr->pmePpCommGpu->getGpuForceStagingPtr()) // buffer received from other GPU
1722 : nullptr; // PME reduction not active on GPU
1724 GpuEventSynchronizer* const pmeSynchronizer =
1725 stepWork.useGpuPmeFReduction
1726 ? (thisRankHasDuty(cr, DUTY_PME) ? pme_gpu_get_f_ready_synchronizer(fr->pmedata)
1727 : // PME force buffer on same GPU
1728 static_cast<GpuEventSynchronizer*>(
1729 fr->pmePpCommGpu->getForcesReadySynchronizer())) // buffer received from other GPU
1730 : nullptr; // PME reduction not active on GPU
1732 gmx::FixedCapacityVector<GpuEventSynchronizer*, 3> dependencyList;
1734 if (stepWork.useGpuPmeFReduction)
1736 dependencyList.push_back(pmeSynchronizer);
1739 gmx::ArrayRef<gmx::RVec> forceWithShift = forceOut.forceWithShiftForces().force();
1741 if (stepWork.useGpuFBufferOps)
1743 // Flag to specify whether the CPU force buffer has contributions to
1744 // local atoms. This depends on whether there are CPU-based force tasks
1745 // or when DD is active the halo exchange has resulted in contributions
1746 // from the non-local part.
1747 const bool haveLocalForceContribInCpuBuffer =
1748 (domainWork.haveCpuLocalForceWork || havePPDomainDecomposition(cr));
1750 // TODO: move these steps as early as possible:
1751 // - CPU f H2D should be as soon as all CPU-side forces are done
1752 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
1753 // before the next CPU task that consumes the forces: vsite spread or update)
1754 // - copy is not perfomed if GPU force halo exchange is active, because it would overwrite the result
1755 // of the halo exchange. In that case the copy is instead performed above, before the exchange.
1756 // These should be unified.
1757 if (haveLocalForceContribInCpuBuffer && !useGpuForcesHaloExchange)
1759 // Note: AtomLocality::All is used for the non-DD case because, as in this
1760 // case copyForcesToGpu() uses a separate stream, it allows overlap of
1761 // CPU force H2D with GPU force tasks on all streams including those in the
1762 // local stream which would otherwise be implicit dependencies for the
1763 // transfer and would not overlap.
1764 auto locality = havePPDomainDecomposition(cr) ? AtomLocality::Local : AtomLocality::All;
1766 stateGpu->copyForcesToGpu(forceWithShift, locality);
1767 dependencyList.push_back(
1768 stateGpu->getForcesReadyOnDeviceEvent(locality, stepWork.useGpuFBufferOps));
1770 if (useGpuForcesHaloExchange)
1772 dependencyList.push_back(cr->dd->gpuHaloExchange[0]->getForcesReadyOnDeviceEvent());
1774 nbv->atomdata_add_nbat_f_to_f_gpu(AtomLocality::Local, stateGpu->getForces(), pmeForcePtr,
1775 dependencyList, stepWork.useGpuPmeFReduction,
1776 haveLocalForceContribInCpuBuffer);
1777 // Copy forces to host if they are needed for update or if virtual sites are enabled.
1778 // If there are vsites, we need to copy forces every step to spread vsite forces on host.
1779 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1780 // copy call done in sim_utils(...) for the output.
1781 // NOTE: If there are virtual sites, the forces are modified on host after this D2H copy. Hence,
1782 // they should not be copied in do_md(...) for the output.
1783 if (!simulationWork.useGpuUpdate || vsite)
1785 stateGpu->copyForcesFromGpu(forceWithShift, AtomLocality::Local);
1786 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1791 nbv->atomdata_add_nbat_f_to_f(AtomLocality::Local, forceWithShift);
1795 launchGpuEndOfStepTasks(nbv, fr->gpuBonded, fr->pmedata, enerd, *runScheduleWork,
1796 useGpuPmeOnThisRank, step, wcycle);
1798 if (DOMAINDECOMP(cr))
1800 dd_force_flop_stop(cr->dd, nrnb);
1803 if (stepWork.computeForces)
1805 rvec* f = as_rvec_array(forceOut.forceWithShiftForces().force().data());
1807 /* If we have NoVirSum forces, but we do not calculate the virial,
1808 * we sum fr->f_novirsum=forceOut.f later.
1810 if (vsite && !(fr->haveDirectVirialContributions && !stepWork.computeVirial))
1812 rvec* fshift = as_rvec_array(forceOut.forceWithShiftForces().shiftForces().data());
1813 spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), f, fshift, FALSE,
1814 nullptr, nrnb, top->idef, fr->pbcType, fr->bMolPBC, graph, box, cr, wcycle);
1817 if (stepWork.computeVirial)
1819 /* Calculation of the virial must be done after vsites! */
1820 calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()),
1821 forceOut.forceWithShiftForces(), vir_force, graph, box, nrnb, fr,
1826 // TODO refactor this and unify with above GPU PME-PP / GPU update path call to the same function
1827 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME) && !simulationWork.useGpuPmePpCommunication
1828 && !simulationWork.useGpuUpdate)
1830 /* In case of node-splitting, the PP nodes receive the long-range
1831 * forces, virial and energy from the PME nodes here.
1833 pme_receive_force_ener(fr, cr, &forceOut.forceWithVirial(), enerd,
1834 simulationWork.useGpuPmePpCommunication, false, wcycle);
1837 if (stepWork.computeForces)
1839 post_process_forces(cr, step, nrnb, wcycle, top, box, as_rvec_array(x.unpaddedArrayRef().data()),
1840 &forceOut, vir_force, mdatoms, graph, fr, vsite, stepWork);
1843 if (stepWork.computeEnergy)
1845 /* Sum the potential energy terms from group contributions */
1846 sum_epot(&(enerd->grpp), enerd->term);
1848 if (!EI_TPI(inputrec->eI))
1850 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1854 /* In case we don't have constraints and are using GPUs, the next balancing
1855 * region starts here.
1856 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1857 * virial calculation and COM pulling, is not thus not included in
1858 * the balance timing, which is ok as most tasks do communication.
1860 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);