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/nonbonded/nb_free_energy.h"
58 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
59 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
60 #include "gromacs/gpu_utils/gpu_utils.h"
61 #include "gromacs/imd/imd.h"
62 #include "gromacs/listed_forces/disre.h"
63 #include "gromacs/listed_forces/gpubonded.h"
64 #include "gromacs/listed_forces/listed_forces.h"
65 #include "gromacs/listed_forces/manage_threading.h"
66 #include "gromacs/listed_forces/orires.h"
67 #include "gromacs/math/arrayrefwithpadding.h"
68 #include "gromacs/math/functions.h"
69 #include "gromacs/math/units.h"
70 #include "gromacs/math/vec.h"
71 #include "gromacs/math/vecdump.h"
72 #include "gromacs/mdlib/calcmu.h"
73 #include "gromacs/mdlib/calcvir.h"
74 #include "gromacs/mdlib/constr.h"
75 #include "gromacs/mdlib/enerdata_utils.h"
76 #include "gromacs/mdlib/force.h"
77 #include "gromacs/mdlib/forcerec.h"
78 #include "gromacs/mdlib/gmx_omp_nthreads.h"
79 #include "gromacs/mdlib/ppforceworkload.h"
80 #include "gromacs/mdlib/qmmm.h"
81 #include "gromacs/mdlib/update.h"
82 #include "gromacs/mdtypes/commrec.h"
83 #include "gromacs/mdtypes/enerdata.h"
84 #include "gromacs/mdtypes/forceoutput.h"
85 #include "gromacs/mdtypes/iforceprovider.h"
86 #include "gromacs/mdtypes/inputrec.h"
87 #include "gromacs/mdtypes/md_enums.h"
88 #include "gromacs/mdtypes/state.h"
89 #include "gromacs/nbnxm/atomdata.h"
90 #include "gromacs/nbnxm/gpu_data_mgmt.h"
91 #include "gromacs/nbnxm/nbnxm.h"
92 #include "gromacs/pbcutil/ishift.h"
93 #include "gromacs/pbcutil/mshift.h"
94 #include "gromacs/pbcutil/pbc.h"
95 #include "gromacs/pulling/pull.h"
96 #include "gromacs/pulling/pull_rotation.h"
97 #include "gromacs/timing/cyclecounter.h"
98 #include "gromacs/timing/gpu_timing.h"
99 #include "gromacs/timing/wallcycle.h"
100 #include "gromacs/timing/wallcyclereporting.h"
101 #include "gromacs/timing/walltime_accounting.h"
102 #include "gromacs/topology/topology.h"
103 #include "gromacs/utility/arrayref.h"
104 #include "gromacs/utility/basedefinitions.h"
105 #include "gromacs/utility/cstringutil.h"
106 #include "gromacs/utility/exceptions.h"
107 #include "gromacs/utility/fatalerror.h"
108 #include "gromacs/utility/gmxassert.h"
109 #include "gromacs/utility/gmxmpi.h"
110 #include "gromacs/utility/logger.h"
111 #include "gromacs/utility/smalloc.h"
112 #include "gromacs/utility/strconvert.h"
113 #include "gromacs/utility/sysinfo.h"
115 using gmx::ForceOutputs;
117 // TODO: this environment variable allows us to verify before release
118 // that on less common architectures the total cost of polling is not larger than
119 // a blocking wait (so polling does not introduce overhead when the static
120 // PME-first ordering would suffice).
121 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
123 // environment variable to enable GPU buffer ops, to allow incremental and optional
124 // introduction of this functionality.
125 // TODO eventially tie this in with other existing GPU flags.
126 static const bool c_enableGpuBufOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
128 static void sum_forces(rvec f[], gmx::ArrayRef<const gmx::RVec> forceToAdd)
130 const int end = forceToAdd.size();
132 int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
133 #pragma omp parallel for num_threads(nt) schedule(static)
134 for (int i = 0; i < end; i++)
136 rvec_inc(f[i], forceToAdd[i]);
140 static void calc_virial(int start, int homenr, const rvec x[],
141 const gmx::ForceWithShiftForces &forceWithShiftForces,
142 tensor vir_part, const t_graph *graph, const matrix box,
143 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
145 /* The short-range virial from surrounding boxes */
146 const rvec *fshift = as_rvec_array(forceWithShiftForces.shiftForces().data());
147 calc_vir(SHIFTS, fr->shift_vec, fshift, vir_part, ePBC == epbcSCREW, box);
148 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
150 /* Calculate partial virial, for local atoms only, based on short range.
151 * Total virial is computed in global_stat, called from do_md
153 const rvec *f = as_rvec_array(forceWithShiftForces.force().data());
154 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
155 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
159 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
163 static void pull_potential_wrapper(const t_commrec *cr,
164 const t_inputrec *ir,
165 const matrix box, gmx::ArrayRef<const gmx::RVec> x,
166 gmx::ForceWithVirial *force,
167 const t_mdatoms *mdatoms,
168 gmx_enerdata_t *enerd,
172 gmx_wallcycle_t wcycle)
177 /* Calculate the center of mass forces, this requires communication,
178 * which is why pull_potential is called close to other communication.
180 wallcycle_start(wcycle, ewcPULLPOT);
181 set_pbc(&pbc, ir->ePBC, box);
183 enerd->term[F_COM_PULL] +=
184 pull_potential(pull_work, mdatoms, &pbc,
185 cr, t, lambda[efptRESTRAINT], as_rvec_array(x.data()), force, &dvdl);
186 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
187 wallcycle_stop(wcycle, ewcPULLPOT);
190 static void pme_receive_force_ener(const t_commrec *cr,
191 gmx::ForceWithVirial *forceWithVirial,
192 gmx_enerdata_t *enerd,
193 gmx_wallcycle_t wcycle)
195 real e_q, e_lj, dvdl_q, dvdl_lj;
196 float cycles_ppdpme, cycles_seppme;
198 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
199 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
201 /* In case of node-splitting, the PP nodes receive the long-range
202 * forces, virial and energy from the PME nodes here.
204 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
207 gmx_pme_receive_f(cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
209 enerd->term[F_COUL_RECIP] += e_q;
210 enerd->term[F_LJ_RECIP] += e_lj;
211 enerd->dvdl_lin[efptCOUL] += dvdl_q;
212 enerd->dvdl_lin[efptVDW] += dvdl_lj;
216 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
218 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
221 static void print_large_forces(FILE *fp,
229 real force2Tolerance = gmx::square(forceTolerance);
230 gmx::index numNonFinite = 0;
231 for (int i = 0; i < md->homenr; i++)
233 real force2 = norm2(f[i]);
234 bool nonFinite = !std::isfinite(force2);
235 if (force2 >= force2Tolerance || nonFinite)
237 fprintf(fp, "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
239 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
246 if (numNonFinite > 0)
248 /* Note that with MPI this fatal call on one rank might interrupt
249 * the printing on other ranks. But we can only avoid that with
250 * an expensive MPI barrier that we would need at each step.
252 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
256 static void post_process_forces(const t_commrec *cr,
259 gmx_wallcycle_t wcycle,
260 const gmx_localtop_t *top,
263 ForceOutputs *forceOutputs,
265 const t_mdatoms *mdatoms,
266 const t_graph *graph,
267 const t_forcerec *fr,
268 const gmx_vsite_t *vsite,
269 const gmx::ForceFlags &forceFlags)
271 rvec *f = as_rvec_array(forceOutputs->forceWithShiftForces().force().data());
273 if (fr->haveDirectVirialContributions)
275 auto &forceWithVirial = forceOutputs->forceWithVirial();
276 rvec *fDirectVir = as_rvec_array(forceWithVirial.force_.data());
280 /* Spread the mesh force on virtual sites to the other particles...
281 * This is parallellized. MPI communication is performed
282 * if the constructing atoms aren't local.
284 matrix virial = { { 0 } };
285 spread_vsite_f(vsite, x, fDirectVir, nullptr,
286 forceFlags.computeVirial, virial,
288 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
289 forceWithVirial.addVirialContribution(virial);
292 if (forceFlags.computeVirial)
294 /* Now add the forces, this is local */
295 sum_forces(f, forceWithVirial.force_);
297 /* Add the direct virial contributions */
298 GMX_ASSERT(forceWithVirial.computeVirial_, "forceWithVirial should request virial computation when we request the virial");
299 m_add(vir_force, forceWithVirial.getVirial(), vir_force);
303 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
308 if (fr->print_force >= 0)
310 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
314 static void do_nb_verlet(t_forcerec *fr,
315 const interaction_const_t *ic,
316 gmx_enerdata_t *enerd,
317 const gmx::ForceFlags &forceFlags,
318 const Nbnxm::InteractionLocality ilocality,
322 gmx_wallcycle_t wcycle)
324 if (!forceFlags.computeNonbondedForces)
326 /* skip non-bonded calculation */
330 nonbonded_verlet_t *nbv = fr->nbv.get();
332 /* GPU kernel launch overhead is already timed separately */
333 if (fr->cutoff_scheme != ecutsVERLET)
335 gmx_incons("Invalid cut-off scheme passed!");
340 /* When dynamic pair-list pruning is requested, we need to prune
341 * at nstlistPrune steps.
343 if (nbv->isDynamicPruningStepCpu(step))
345 /* Prune the pair-list beyond fr->ic->rlistPrune using
346 * the current coordinates of the atoms.
348 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
349 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
350 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
354 nbv->dispatchNonbondedKernel(ilocality, *ic, forceFlags, clearF, *fr, enerd, nrnb);
357 static inline void clear_rvecs_omp(int n, rvec v[])
359 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, n);
361 /* Note that we would like to avoid this conditional by putting it
362 * into the omp pragma instead, but then we still take the full
363 * omp parallel for overhead (at least with gcc5).
367 for (int i = 0; i < n; i++)
374 #pragma omp parallel for num_threads(nth) schedule(static)
375 for (int i = 0; i < n; i++)
382 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
384 * \param groupOptions Group options, containing T-coupling options
386 static real averageKineticEnergyEstimate(const t_grpopts &groupOptions)
388 real nrdfCoupled = 0;
389 real nrdfUncoupled = 0;
390 real kineticEnergy = 0;
391 for (int g = 0; g < groupOptions.ngtc; g++)
393 if (groupOptions.tau_t[g] >= 0)
395 nrdfCoupled += groupOptions.nrdf[g];
396 kineticEnergy += groupOptions.nrdf[g]*0.5*groupOptions.ref_t[g]*BOLTZ;
400 nrdfUncoupled += groupOptions.nrdf[g];
404 /* This conditional with > also catches nrdf=0 */
405 if (nrdfCoupled > nrdfUncoupled)
407 return kineticEnergy*(nrdfCoupled + nrdfUncoupled)/nrdfCoupled;
415 /*! \brief This routine checks that the potential energy is finite.
417 * Always checks that the potential energy is finite. If step equals
418 * inputrec.init_step also checks that the magnitude of the potential energy
419 * is reasonable. Terminates with a fatal error when a check fails.
420 * Note that passing this check does not guarantee finite forces,
421 * since those use slightly different arithmetics. But in most cases
422 * there is just a narrow coordinate range where forces are not finite
423 * and energies are finite.
425 * \param[in] step The step number, used for checking and printing
426 * \param[in] enerd The energy data; the non-bonded group energies need to be added to enerd.term[F_EPOT] before calling this routine
427 * \param[in] inputrec The input record
429 static void checkPotentialEnergyValidity(int64_t step,
430 const gmx_enerdata_t &enerd,
431 const t_inputrec &inputrec)
433 /* Threshold valid for comparing absolute potential energy against
434 * the kinetic energy. Normally one should not consider absolute
435 * potential energy values, but with a factor of one million
436 * we should never get false positives.
438 constexpr real c_thresholdFactor = 1e6;
440 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
441 real averageKineticEnergy = 0;
442 /* We only check for large potential energy at the initial step,
443 * because that is by far the most likely step for this too occur
444 * and because computing the average kinetic energy is not free.
445 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
446 * before they become NaN.
448 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
450 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
453 if (energyIsNotFinite || (averageKineticEnergy > 0 &&
454 enerd.term[F_EPOT] > c_thresholdFactor*averageKineticEnergy))
456 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.",
459 energyIsNotFinite ? "not finite" : "extremely high",
461 enerd.term[F_COUL_SR],
462 energyIsNotFinite ? "non-finite" : "very high",
463 energyIsNotFinite ? " or Nan" : "");
467 /*! \brief Return true if there are special forces computed this step.
469 * The conditionals exactly correspond to those in computeSpecialForces().
472 haveSpecialForces(const t_inputrec *inputrec,
473 ForceProviders *forceProviders,
474 const pull_t *pull_work,
475 const bool computeForces,
480 ((computeForces && forceProviders->hasForceProvider()) || // forceProviders
481 (inputrec->bPull && pull_have_potential(pull_work)) || // pull
482 inputrec->bRot || // enforced rotation
483 (ed != nullptr) || // flooding
484 (inputrec->bIMD && computeForces)); // IMD
487 /*! \brief Compute forces and/or energies for special algorithms
489 * The intention is to collect all calls to algorithms that compute
490 * forces on local atoms only and that do not contribute to the local
491 * virial sum (but add their virial contribution separately).
492 * Eventually these should likely all become ForceProviders.
493 * Within this function the intention is to have algorithms that do
494 * global communication at the end, so global barriers within the MD loop
495 * are as close together as possible.
497 * \param[in] fplog The log file
498 * \param[in] cr The communication record
499 * \param[in] inputrec The input record
500 * \param[in] awh The Awh module (nullptr if none in use).
501 * \param[in] enforcedRotation Enforced rotation module.
502 * \param[in] imdSession The IMD session
503 * \param[in] pull_work The pull work structure.
504 * \param[in] step The current MD step
505 * \param[in] t The current time
506 * \param[in,out] wcycle Wallcycle accounting struct
507 * \param[in,out] forceProviders Pointer to a list of force providers
508 * \param[in] box The unit cell
509 * \param[in] x The coordinates
510 * \param[in] mdatoms Per atom properties
511 * \param[in] lambda Array of free-energy lambda values
512 * \param[in] forceFlags Force schedule flags
513 * \param[in,out] forceWithVirial Force and virial buffers
514 * \param[in,out] enerd Energy buffer
515 * \param[in,out] ed Essential dynamics pointer
516 * \param[in] didNeighborSearch Tells if we did neighbor searching this step, used for ED sampling
518 * \todo Remove didNeighborSearch, which is used incorrectly.
519 * \todo Convert all other algorithms called here to ForceProviders.
522 computeSpecialForces(FILE *fplog,
524 const t_inputrec *inputrec,
526 gmx_enfrot *enforcedRotation,
527 gmx::ImdSession *imdSession,
531 gmx_wallcycle_t wcycle,
532 ForceProviders *forceProviders,
534 gmx::ArrayRef<const gmx::RVec> x,
535 const t_mdatoms *mdatoms,
537 const gmx::ForceFlags &forceFlags,
538 gmx::ForceWithVirial *forceWithVirial,
539 gmx_enerdata_t *enerd,
541 bool didNeighborSearch)
543 /* NOTE: Currently all ForceProviders only provide forces.
544 * When they also provide energies, remove this conditional.
546 if (forceFlags.computeForces)
548 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
549 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
551 /* Collect forces from modules */
552 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
555 if (inputrec->bPull && pull_have_potential(pull_work))
557 pull_potential_wrapper(cr, inputrec, box, x,
559 mdatoms, enerd, pull_work, lambda, t,
564 enerd->term[F_COM_PULL] +=
565 awh->applyBiasForcesAndUpdateBias(inputrec->ePBC, *mdatoms, box,
567 t, step, wcycle, fplog);
571 rvec *f = as_rvec_array(forceWithVirial->force_.data());
573 /* Add the forces from enforced rotation potentials (if any) */
576 wallcycle_start(wcycle, ewcROTadd);
577 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
578 wallcycle_stop(wcycle, ewcROTadd);
583 /* Note that since init_edsam() is called after the initialization
584 * of forcerec, edsam doesn't request the noVirSum force buffer.
585 * Thus if no other algorithm (e.g. PME) requires it, the forces
586 * here will contribute to the virial.
588 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, didNeighborSearch);
591 /* Add forces from interactive molecular dynamics (IMD), if any */
592 if (inputrec->bIMD && forceFlags.computeForces)
594 imdSession->applyForces(f);
598 /*! \brief Launch the prepare_step and spread stages of PME GPU.
600 * \param[in] pmedata The PME structure
601 * \param[in] box The box matrix
602 * \param[in] x Coordinate array
603 * \param[in] forceFlags Force schedule flags
604 * \param[in] pmeFlags PME flags
605 * \param[in] wcycle The wallcycle structure
607 static inline void launchPmeGpuSpread(gmx_pme_t *pmedata,
610 const gmx::ForceFlags &forceFlags,
612 gmx_wallcycle_t wcycle)
614 pme_gpu_prepare_computation(pmedata, forceFlags.haveDynamicBox, box, wcycle, pmeFlags);
615 pme_gpu_copy_coordinates_to_gpu(pmedata, x, wcycle);
616 pme_gpu_launch_spread(pmedata, wcycle);
619 /*! \brief Launch the FFT and gather stages of PME GPU
621 * This function only implements setting the output forces (no accumulation).
623 * \param[in] pmedata The PME structure
624 * \param[in] wcycle The wallcycle structure
625 * \param[in] useGpuFPmeReduction Whether forces will be reduced on GPU
627 static void launchPmeGpuFftAndGather(gmx_pme_t *pmedata,
628 gmx_wallcycle_t wcycle,
629 bool useGpuFPmeReduction)
631 pme_gpu_launch_complex_transforms(pmedata, wcycle);
632 pme_gpu_launch_gather(pmedata, wcycle, PmeForceOutputHandling::Set, useGpuFPmeReduction);
636 * Polling wait for either of the PME or nonbonded GPU tasks.
638 * Instead of a static order in waiting for GPU tasks, this function
639 * polls checking which of the two tasks completes first, and does the
640 * associated force buffer reduction overlapped with the other task.
641 * By doing that, unlike static scheduling order, it can always overlap
642 * one of the reductions, regardless of the GPU task completion order.
644 * \param[in] nbv Nonbonded verlet structure
645 * \param[in,out] pmedata PME module data
646 * \param[in,out] forceOutputs Output buffer for the forces and virial
647 * \param[in,out] enerd Energy data structure results are reduced into
648 * \param[in] forceFlags Force schedule flags
649 * \param[in] pmeFlags PME flags
650 * \param[in] wcycle The wallcycle structure
652 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t *nbv,
654 gmx::ForceOutputs *forceOutputs,
655 gmx_enerdata_t *enerd,
656 const gmx::ForceFlags &forceFlags,
658 gmx_wallcycle_t wcycle)
660 bool isPmeGpuDone = false;
661 bool isNbGpuDone = false;
665 gmx::ForceWithShiftForces &forceWithShiftForces = forceOutputs->forceWithShiftForces();
666 gmx::ForceWithVirial &forceWithVirial = forceOutputs->forceWithVirial();
668 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
670 while (!isPmeGpuDone || !isNbGpuDone)
674 GpuTaskCompletion completionType = (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
675 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, pmeFlags, wcycle, &forceWithVirial, enerd, completionType);
680 GpuTaskCompletion completionType = (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
681 isNbGpuDone = Nbnxm::gpu_try_finish_task(nbv->gpu_nbv,
683 Nbnxm::AtomLocality::Local,
684 enerd->grpp.ener[egLJSR].data(),
685 enerd->grpp.ener[egCOULSR].data(),
686 forceWithShiftForces.shiftForces(), completionType, wcycle);
690 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
691 forceWithShiftForces.force());
697 /*! \brief Set up the different force buffers; also does clearing.
699 * \param[in] fr force record pointer
700 * \param[in] pull_work The pull work object.
701 * \param[in] inputrec input record
702 * \param[in] force force array
703 * \param[in] forceFlags Force schedule flags
704 * \param[out] wcycle wallcycle recording structure
706 * \returns Cleared force output structure
709 setupForceOutputs(t_forcerec *fr,
711 const t_inputrec &inputrec,
712 gmx::ArrayRefWithPadding<gmx::RVec> force,
713 const gmx::ForceFlags &forceFlags,
714 gmx_wallcycle_t wcycle)
716 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
718 /* NOTE: We assume fr->shiftForces is all zeros here */
719 gmx::ForceWithShiftForces forceWithShiftForces(force, forceFlags.computeVirial, fr->shiftForces);
721 if (forceFlags.computeForces)
723 /* Clear the short- and long-range forces */
724 clear_rvecs_omp(fr->natoms_force_constr,
725 as_rvec_array(forceWithShiftForces.force().data()));
728 /* If we need to compute the virial, we might need a separate
729 * force buffer for algorithms for which the virial is calculated
730 * directly, such as PME. Otherwise, forceWithVirial uses the
731 * the same force (f in legacy calls) buffer as other algorithms.
733 const bool useSeparateForceWithVirialBuffer = (forceFlags.computeForces &&
734 (forceFlags.computeVirial && fr->haveDirectVirialContributions));
735 /* forceWithVirial uses the local atom range only */
736 gmx::ForceWithVirial forceWithVirial(useSeparateForceWithVirialBuffer ?
737 fr->forceBufferForDirectVirialContributions : force.unpaddedArrayRef(),
738 forceFlags.computeVirial);
740 if (useSeparateForceWithVirialBuffer)
742 /* TODO: update comment
743 * We only compute forces on local atoms. Note that vsites can
744 * spread to non-local atoms, but that part of the buffer is
745 * cleared separately in the vsite spreading code.
747 clear_rvecs_omp(forceWithVirial.force_.size(), as_rvec_array(forceWithVirial.force_.data()));
750 if (inputrec.bPull && pull_have_constraint(pull_work))
752 clear_pull_forces(pull_work);
755 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
757 return ForceOutputs(forceWithShiftForces, forceWithVirial);
761 /*! \brief Set up flags that indicate what type of work is there to compute.
763 * Currently we only update it at search steps,
764 * but some properties may change more frequently (e.g. virial/non-virial step),
765 * so when including those either the frequency of update (per-step) or the scope
766 * of a flag will change (i.e. a set of flags for nstlist steps).
770 setupForceWorkload(gmx::PpForceWorkload *forceWork,
771 const t_inputrec *inputrec,
772 const t_forcerec *fr,
773 const pull_t *pull_work,
777 const gmx::ForceFlags &forceFlags
780 forceWork->haveSpecialForces = haveSpecialForces(inputrec, fr->forceProviders, pull_work, forceFlags.computeForces, ed);
781 forceWork->haveCpuBondedWork = haveCpuBondeds(*fr);
782 forceWork->haveGpuBondedWork = ((fr->gpuBonded != nullptr) && fr->gpuBonded->haveInteractions());
783 forceWork->haveRestraintsWork = havePositionRestraints(idef, *fcd);
784 forceWork->haveCpuListedForceWork = haveCpuListedForces(*fr, idef, *fcd);
787 /*! \brief Set up force flag stuct from the force bitmask.
789 * \param[out] flags Force schedule flags
790 * \param[in] legacyFlags Force bitmask flags used to construct the new flags
791 * \param[in] isNonbondedOn Global override, if false forces to turn off all nonbonded calculation.
794 setupForceFlags(gmx::ForceFlags *flags,
795 const int legacyFlags,
796 const bool isNonbondedOn)
798 flags->stateChanged = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
799 flags->haveDynamicBox = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
800 flags->doNeighborSearch = ((legacyFlags & GMX_FORCE_NS) != 0);
801 flags->computeVirial = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
802 flags->computeEnergy = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
803 flags->computeForces = ((legacyFlags & GMX_FORCE_FORCES) != 0);
804 flags->computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
805 flags->computeNonbondedForces = ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && isNonbondedOn;
806 flags->computeDhdl = ((legacyFlags & GMX_FORCE_DHDL) != 0);
810 /* \brief Launch end-of-step GPU tasks: buffer clearing and rolling pruning.
812 * TODO: eliminate the \p useGpuNonbonded and \p useGpuNonbonded when these are
813 * incorporated in PpForceWorkload.
816 launchGpuEndOfStepTasks(nonbonded_verlet_t *nbv,
817 gmx::GpuBonded *gpuBonded,
819 gmx_enerdata_t *enerd,
820 const gmx::MdScheduleWorkload &mdScheduleWork,
821 bool useGpuNonbonded,
824 gmx_wallcycle_t wcycle)
828 /* Launch pruning before buffer clearing because the API overhead of the
829 * clear kernel launches can leave the GPU idle while it could be running
832 if (nbv->isDynamicPruningStepGpu(step))
834 nbv->dispatchPruneKernelGpu(step);
837 /* now clear the GPU outputs while we finish the step on the CPU */
838 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
839 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
840 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, mdScheduleWork.forceFlags.computeVirial);
841 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
842 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
847 pme_gpu_reinit_computation(pmedata, wcycle);
850 if (mdScheduleWork.forceWork.haveGpuBondedWork && mdScheduleWork.forceFlags.computeEnergy)
852 // in principle this should be included in the DD balancing region,
853 // but generally it is infrequent so we'll omit it for the sake of
855 gpuBonded->waitAccumulateEnergyTerms(enerd);
857 gpuBonded->clearEnergies();
862 void do_force(FILE *fplog,
864 const gmx_multisim_t *ms,
865 const t_inputrec *inputrec,
867 gmx_enfrot *enforcedRotation,
868 gmx::ImdSession *imdSession,
872 gmx_wallcycle_t wcycle,
873 const gmx_localtop_t *top,
875 gmx::ArrayRefWithPadding<gmx::RVec> x,
877 gmx::ArrayRefWithPadding<gmx::RVec> force,
879 const t_mdatoms *mdatoms,
880 gmx_enerdata_t *enerd,
882 gmx::ArrayRef<real> lambda,
885 gmx::MdScheduleWorkload *mdScheduleWork,
886 const gmx_vsite_t *vsite,
891 const DDBalanceRegionHandler &ddBalanceRegionHandler)
895 gmx_bool bFillGrid, bCalcCGCM;
896 gmx_bool bUseGPU, bUseOrEmulGPU;
897 nonbonded_verlet_t *nbv = fr->nbv.get();
898 interaction_const_t *ic = fr->ic;
900 // TODO remove the code below when the legacy flags are not in use anymore
901 /* modify force flag if not doing nonbonded */
904 legacyFlags &= ~GMX_FORCE_NONBONDED;
906 setupForceFlags(&mdScheduleWork->forceFlags, legacyFlags, fr->bNonbonded);
908 const gmx::ForceFlags &forceFlags = mdScheduleWork->forceFlags;
910 bFillGrid = (forceFlags.doNeighborSearch && forceFlags.stateChanged);
911 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
912 bUseGPU = fr->nbv->useGpu();
913 bUseOrEmulGPU = bUseGPU || fr->nbv->emulateGpu();
915 const auto pmeRunMode = fr->pmedata ? pme_run_mode(fr->pmedata) : PmeRunMode::CPU;
916 // TODO slim this conditional down - inputrec and duty checks should mean the same in proper code!
917 const bool useGpuPme = EEL_PME(fr->ic->eeltype) && thisRankHasDuty(cr, DUTY_PME) &&
918 ((pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed));
919 const int pmeFlags = GMX_PME_SPREAD | GMX_PME_SOLVE |
920 (forceFlags.computeVirial ? GMX_PME_CALC_ENER_VIR : 0) |
921 (forceFlags.computeEnergy ? GMX_PME_CALC_ENER_VIR : 0) |
922 (forceFlags.computeForces ? GMX_PME_CALC_F : 0);
924 // Switches on whether to use GPU for position and force buffer operations
925 // TODO consider all possible combinations of triggers, and how to combine optimally in each case.
926 const BufferOpsUseGpu useGpuXBufOps = (c_enableGpuBufOps && bUseGPU && (GMX_GPU == GMX_GPU_CUDA)) ?
927 BufferOpsUseGpu::True : BufferOpsUseGpu::False;;
928 // GPU Force buffer ops are disabled on virial steps, because the virial calc is not yet ported to GPU
929 const BufferOpsUseGpu useGpuFBufOps = (c_enableGpuBufOps && bUseGPU && (GMX_GPU == GMX_GPU_CUDA))
930 && !(forceFlags.computeVirial || forceFlags.computeEnergy) ?
931 BufferOpsUseGpu::True : BufferOpsUseGpu::False;
932 // TODO: move / add this flag to the internal PME GPU data structures
933 const bool useGpuFPmeReduction = (useGpuFBufOps == BufferOpsUseGpu::True) &&
934 thisRankHasDuty(cr, DUTY_PME) && useGpuPme; // only supported if this rank is perfoming PME on the GPU
936 /* At a search step we need to start the first balancing region
937 * somewhere early inside the step after communication during domain
938 * decomposition (and not during the previous step as usual).
940 if (forceFlags.doNeighborSearch)
942 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
946 const int homenr = mdatoms->homenr;
948 clear_mat(vir_force);
950 if (forceFlags.stateChanged)
952 if (inputrecNeedMutot(inputrec))
954 /* Calculate total (local) dipole moment in a temporary common array.
955 * This makes it possible to sum them over nodes faster.
957 calc_mu(start, homenr,
958 x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
963 if (fr->ePBC != epbcNONE)
965 /* Compute shift vectors every step,
966 * because of pressure coupling or box deformation!
968 if (forceFlags.haveDynamicBox && forceFlags.stateChanged)
970 calc_shifts(box, fr->shift_vec);
975 put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr), gmx_omp_nthreads_get(emntDefault));
976 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
978 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
980 unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
984 nbnxn_atomdata_copy_shiftvec(forceFlags.haveDynamicBox,
985 fr->shift_vec, nbv->nbat.get());
988 if (!thisRankHasDuty(cr, DUTY_PME))
990 /* Send particle coordinates to the pme nodes.
991 * Since this is only implemented for domain decomposition
992 * and domain decomposition does not use the graph,
993 * we do not need to worry about shifting.
995 gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
996 lambda[efptCOUL], lambda[efptVDW],
997 (forceFlags.computeVirial || forceFlags.computeEnergy),
1000 #endif /* GMX_MPI */
1004 launchPmeGpuSpread(fr->pmedata, box, as_rvec_array(x.unpaddedArrayRef().data()), forceFlags, pmeFlags, wcycle);
1007 /* do gridding for pair search */
1008 if (forceFlags.doNeighborSearch)
1010 if (graph && forceFlags.stateChanged)
1012 /* Calculate intramolecular shift vectors to make molecules whole */
1013 mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
1017 // - vzero is constant, do we need to pass it?
1018 // - box_diag should be passed directly to nbnxn_put_on_grid
1024 box_diag[XX] = box[XX][XX];
1025 box_diag[YY] = box[YY][YY];
1026 box_diag[ZZ] = box[ZZ][ZZ];
1028 wallcycle_start(wcycle, ewcNS);
1029 if (!DOMAINDECOMP(cr))
1031 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
1032 nbnxn_put_on_grid(nbv, box,
1034 nullptr, 0, mdatoms->homenr, -1,
1035 fr->cginfo, x.unpaddedArrayRef(),
1037 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
1041 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
1042 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd),
1043 fr->cginfo, x.unpaddedArrayRef());
1044 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
1047 nbv->setAtomProperties(*mdatoms, fr->cginfo);
1049 wallcycle_stop(wcycle, ewcNS);
1051 /* initialize the GPU nbnxm atom data and bonded data structures */
1054 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1056 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1057 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
1058 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1062 /* Now we put all atoms on the grid, we can assign bonded
1063 * interactions to the GPU, where the grid order is
1064 * needed. Also the xq, f and fshift device buffers have
1065 * been reallocated if needed, so the bonded code can
1066 * learn about them. */
1067 // TODO the xq, f, and fshift buffers are now shared
1068 // resources, so they should be maintained by a
1069 // higher-level object than the nb module.
1070 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1072 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1073 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1074 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1076 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1080 // Call it per-step as force-flags can change.
1081 // Need to run after the GPU-offload bonded interaction lists
1082 // are set up to be able to determine whether there is bonded work.
1083 setupForceWorkload(&mdScheduleWork->forceWork,
1092 const gmx::PpForceWorkload &forceWork = mdScheduleWork->forceWork;
1094 /* do local pair search */
1095 if (forceFlags.doNeighborSearch)
1097 // TODO: fuse this branch with the above forceFlags.doNeighborSearch block
1098 wallcycle_start_nocount(wcycle, ewcNS);
1099 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1100 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1101 nbv->constructPairlist(Nbnxm::InteractionLocality::Local,
1102 &top->excls, step, nrnb);
1104 nbv->setupGpuShortRangeWork(fr->gpuBonded, Nbnxm::InteractionLocality::Local);
1106 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1107 wallcycle_stop(wcycle, ewcNS);
1109 if (useGpuXBufOps == BufferOpsUseGpu::True)
1111 nbv->atomdata_init_copy_x_to_nbat_x_gpu();
1113 // For force buffer ops, we use the below conditon rather than
1114 // useGpuFBufOps to ensure that init is performed even if this
1115 // NS step is also a virial step (on which f buf ops are deactivated).
1116 if (c_enableGpuBufOps && bUseGPU && (GMX_GPU == GMX_GPU_CUDA))
1118 nbv->atomdata_init_add_nbat_f_to_f_gpu();
1123 if (useGpuXBufOps == BufferOpsUseGpu::True)
1125 // The condition here was (pme != nullptr && pme_gpu_get_device_x(fr->pmedata) != nullptr)
1128 nbv->copyCoordinatesToGpu(Nbnxm::AtomLocality::Local, false,
1129 x.unpaddedArrayRef());
1131 nbv->convertCoordinatesGpu(Nbnxm::AtomLocality::Local, false,
1132 useGpuPme ? pme_gpu_get_device_x(fr->pmedata) : nbv->getDeviceCoordinates());
1136 nbv->convertCoordinates(Nbnxm::AtomLocality::Local, false,
1137 x.unpaddedArrayRef());
1143 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1145 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1147 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1148 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1149 if (forceFlags.doNeighborSearch || (useGpuXBufOps == BufferOpsUseGpu::False))
1151 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1152 Nbnxm::AtomLocality::Local);
1154 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1155 // with X buffer ops offloaded to the GPU on all but the search steps
1157 // bonded work not split into separate local and non-local, so with DD
1158 // we can only launch the kernel after non-local coordinates have been received.
1159 if (forceWork.haveGpuBondedWork && !havePPDomainDecomposition(cr))
1161 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1162 fr->gpuBonded->launchKernel(fr, forceFlags, box);
1163 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1166 /* launch local nonbonded work on GPU */
1167 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1168 do_nb_verlet(fr, ic, enerd, forceFlags, Nbnxm::InteractionLocality::Local, enbvClearFNo,
1169 step, nrnb, wcycle);
1170 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1171 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1176 // In PME GPU and mixed mode we launch FFT / gather after the
1177 // X copy/transform to allow overlap as well as after the GPU NB
1178 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1179 // the nonbonded kernel.
1180 launchPmeGpuFftAndGather(fr->pmedata, wcycle, useGpuFPmeReduction);
1183 /* Communicate coordinates and sum dipole if necessary +
1184 do non-local pair search */
1185 if (havePPDomainDecomposition(cr))
1187 if (forceFlags.doNeighborSearch)
1189 // TODO: fuse this branch with the above large forceFlags.doNeighborSearch block
1190 wallcycle_start_nocount(wcycle, ewcNS);
1191 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1192 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1193 nbv->constructPairlist(Nbnxm::InteractionLocality::NonLocal,
1194 &top->excls, step, nrnb);
1196 nbv->setupGpuShortRangeWork(fr->gpuBonded, Nbnxm::InteractionLocality::NonLocal);
1197 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1198 wallcycle_stop(wcycle, ewcNS);
1202 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1204 if (useGpuXBufOps == BufferOpsUseGpu::True)
1206 // The condition here was (pme != nullptr && pme_gpu_get_device_x(fr->pmedata) != nullptr)
1209 nbv->copyCoordinatesToGpu(Nbnxm::AtomLocality::NonLocal, false,
1210 x.unpaddedArrayRef());
1212 nbv->convertCoordinatesGpu(Nbnxm::AtomLocality::NonLocal, false,
1213 useGpuPme ? pme_gpu_get_device_x(fr->pmedata) : nbv->getDeviceCoordinates());
1217 nbv->convertCoordinates(Nbnxm::AtomLocality::NonLocal, false,
1218 x.unpaddedArrayRef());
1225 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1227 if (forceFlags.doNeighborSearch || (useGpuXBufOps == BufferOpsUseGpu::False))
1229 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1230 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1231 Nbnxm::AtomLocality::NonLocal);
1232 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1235 if (forceWork.haveGpuBondedWork)
1237 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1238 fr->gpuBonded->launchKernel(fr, forceFlags, box);
1239 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1242 /* launch non-local nonbonded tasks on GPU */
1243 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1244 do_nb_verlet(fr, ic, enerd, forceFlags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1245 step, nrnb, wcycle);
1246 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1248 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1254 /* launch D2H copy-back F */
1255 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1256 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1258 bool copyBackNbForce = (useGpuFBufOps == BufferOpsUseGpu::False);
1260 if (havePPDomainDecomposition(cr))
1262 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1263 forceFlags, Nbnxm::AtomLocality::NonLocal, copyBackNbForce);
1265 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1266 forceFlags, Nbnxm::AtomLocality::Local, copyBackNbForce);
1267 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1269 if (forceWork.haveGpuBondedWork && forceFlags.computeEnergy)
1271 fr->gpuBonded->launchEnergyTransfer();
1273 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1276 if (forceFlags.stateChanged && inputrecNeedMutot(inputrec))
1280 gmx_sumd(2*DIM, mu, cr);
1282 ddBalanceRegionHandler.reopenRegionCpu();
1285 for (i = 0; i < 2; i++)
1287 for (j = 0; j < DIM; j++)
1289 fr->mu_tot[i][j] = mu[i*DIM + j];
1293 if (fr->efep == efepNO)
1295 copy_rvec(fr->mu_tot[0], mu_tot);
1299 for (j = 0; j < DIM; j++)
1302 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1303 lambda[efptCOUL]*fr->mu_tot[1][j];
1307 /* Reset energies */
1308 reset_enerdata(enerd);
1309 /* Clear the shift forces */
1310 // TODO: This should be linked to the shift force buffer in use, or cleared before use instead
1311 for (gmx::RVec &elem : fr->shiftForces)
1313 elem = { 0.0_real, 0.0_real, 0.0_real };
1316 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1318 wallcycle_start(wcycle, ewcPPDURINGPME);
1319 dd_force_flop_start(cr->dd, nrnb);
1324 wallcycle_start(wcycle, ewcROT);
1325 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, forceFlags.doNeighborSearch);
1326 wallcycle_stop(wcycle, ewcROT);
1329 /* Start the force cycle counter.
1330 * Note that a different counter is used for dynamic load balancing.
1332 wallcycle_start(wcycle, ewcFORCE);
1334 // Set up and clear force outputs.
1335 // We use std::move to keep the compiler happy, it has no effect.
1336 ForceOutputs forceOut = setupForceOutputs(fr, pull_work, *inputrec, std::move(force), forceFlags, wcycle);
1338 /* We calculate the non-bonded forces, when done on the CPU, here.
1339 * We do this before calling do_force_lowlevel, because in that
1340 * function, the listed forces are calculated before PME, which
1341 * does communication. With this order, non-bonded and listed
1342 * force calculation imbalance can be balanced out by the domain
1343 * decomposition load balancing.
1348 do_nb_verlet(fr, ic, enerd, forceFlags, Nbnxm::InteractionLocality::Local, enbvClearFYes,
1349 step, nrnb, wcycle);
1352 if (fr->efep != efepNO)
1354 /* Calculate the local and non-local free energy interactions here.
1355 * Happens here on the CPU both with and without GPU.
1357 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::Local,
1358 fr, as_rvec_array(x.unpaddedArrayRef().data()), &forceOut.forceWithShiftForces(), *mdatoms,
1359 inputrec->fepvals, lambda.data(),
1360 enerd, forceFlags, nrnb);
1362 if (havePPDomainDecomposition(cr))
1364 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::NonLocal,
1365 fr, as_rvec_array(x.unpaddedArrayRef().data()), &forceOut.forceWithShiftForces(), *mdatoms,
1366 inputrec->fepvals, lambda.data(),
1367 enerd, forceFlags, nrnb);
1373 if (havePPDomainDecomposition(cr))
1375 do_nb_verlet(fr, ic, enerd, forceFlags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1376 step, nrnb, wcycle);
1379 /* Add all the non-bonded force to the normal force array.
1380 * This can be split into a local and a non-local part when overlapping
1381 * communication with calculation with domain decomposition.
1383 wallcycle_stop(wcycle, ewcFORCE);
1384 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::All, forceOut.forceWithShiftForces().force());
1385 wallcycle_start_nocount(wcycle, ewcFORCE);
1387 /* If there are multiple fshift output buffers we need to reduce them */
1388 if (forceFlags.computeVirial)
1390 /* This is not in a subcounter because it takes a
1391 negligible and constant-sized amount of time */
1392 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat,
1393 forceOut.forceWithShiftForces().shiftForces());
1397 /* update QMMMrec, if necessary */
1400 update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1403 /* Compute the bonded and non-bonded energies and optionally forces */
1404 do_force_lowlevel(fr, inputrec, &(top->idef),
1405 cr, ms, nrnb, wcycle, mdatoms,
1406 x, hist, &forceOut, enerd, fcd,
1407 box, lambda.data(), graph, fr->mu_tot,
1409 ddBalanceRegionHandler);
1411 wallcycle_stop(wcycle, ewcFORCE);
1413 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1414 imdSession, pull_work, step, t, wcycle,
1415 fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda.data(),
1416 forceFlags, &forceOut.forceWithVirial(), enerd,
1417 ed, forceFlags.doNeighborSearch);
1419 bool useCpuFPmeReduction = thisRankHasDuty(cr, DUTY_PME) && !useGpuFPmeReduction;
1420 bool haveCpuForces = (forceWork.haveSpecialForces || forceWork.haveCpuListedForceWork || useCpuFPmeReduction);
1422 // Will store the amount of cycles spent waiting for the GPU that
1423 // will be later used in the DLB accounting.
1424 float cycles_wait_gpu = 0;
1427 auto &forceWithShiftForces = forceOut.forceWithShiftForces();
1428 rvec *f = as_rvec_array(forceWithShiftForces.force().data());
1430 /* wait for non-local forces (or calculate in emulation mode) */
1431 if (havePPDomainDecomposition(cr))
1435 cycles_wait_gpu += Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1436 forceFlags, Nbnxm::AtomLocality::NonLocal,
1437 enerd->grpp.ener[egLJSR].data(),
1438 enerd->grpp.ener[egCOULSR].data(),
1439 forceWithShiftForces.shiftForces(),
1444 wallcycle_start_nocount(wcycle, ewcFORCE);
1445 do_nb_verlet(fr, ic, enerd, forceFlags, Nbnxm::InteractionLocality::NonLocal, enbvClearFYes,
1446 step, nrnb, wcycle);
1447 wallcycle_stop(wcycle, ewcFORCE);
1450 if (useGpuFBufOps == BufferOpsUseGpu::True && haveCpuForces)
1452 nbv->launch_copy_f_to_gpu(f, Nbnxm::AtomLocality::NonLocal);
1455 // flag to specify if forces should be accumulated in force buffer
1456 // ops. For non-local part, this just depends on whether CPU forces are present.
1457 bool accumulateForce = (useGpuFBufOps == BufferOpsUseGpu::True) && haveCpuForces;
1458 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::NonLocal,
1459 forceWithShiftForces.force(), pme_gpu_get_device_f(fr->pmedata),
1460 pme_gpu_get_f_ready_synchronizer(fr->pmedata),
1461 useGpuFBufOps, useGpuFPmeReduction, accumulateForce);
1462 if (useGpuFBufOps == BufferOpsUseGpu::True)
1464 nbv->launch_copy_f_from_gpu(f, Nbnxm::AtomLocality::NonLocal);
1467 if (fr->nbv->emulateGpu() && forceFlags.computeVirial)
1469 nbnxn_atomdata_add_nbat_fshift_to_fshift(*nbv->nbat,
1470 forceWithShiftForces.shiftForces());
1475 if (havePPDomainDecomposition(cr))
1477 /* We are done with the CPU compute.
1478 * We will now communicate the non-local forces.
1479 * If we use a GPU this will overlap with GPU work, so in that case
1480 * we do not close the DD force balancing region here.
1482 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1484 if (forceFlags.computeForces)
1486 if (useGpuFBufOps == BufferOpsUseGpu::True)
1488 nbv->wait_for_gpu_force_reduction(Nbnxm::AtomLocality::NonLocal);
1490 dd_move_f(cr->dd, &forceOut.forceWithShiftForces(), wcycle);
1494 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1495 // an alternating wait/reduction scheme.
1496 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPme && bUseGPU && !DOMAINDECOMP(cr) &&
1497 (useGpuFBufOps == BufferOpsUseGpu::False));
1498 if (alternateGpuWait)
1500 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &forceOut, enerd,
1501 forceFlags, pmeFlags, wcycle);
1504 if (!alternateGpuWait && useGpuPme)
1506 pme_gpu_wait_and_reduce(fr->pmedata, pmeFlags, wcycle, &forceOut.forceWithVirial(), enerd, useGpuFPmeReduction);
1509 /* Wait for local GPU NB outputs on the non-alternating wait path */
1510 if (!alternateGpuWait && bUseGPU)
1512 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1513 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1514 * but even with a step of 0.1 ms the difference is less than 1%
1517 const float gpuWaitApiOverheadMargin = 2e6F; /* cycles */
1518 const float waitCycles =
1519 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1520 forceFlags, Nbnxm::AtomLocality::Local,
1521 enerd->grpp.ener[egLJSR].data(),
1522 enerd->grpp.ener[egCOULSR].data(),
1523 forceOut.forceWithShiftForces().shiftForces(),
1526 if (ddBalanceRegionHandler.useBalancingRegion())
1528 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1529 if (forceFlags.computeForces && waitCycles <= gpuWaitApiOverheadMargin)
1531 /* We measured few cycles, it could be that the kernel
1532 * and transfer finished earlier and there was no actual
1533 * wait time, only API call overhead.
1534 * Then the actual time could be anywhere between 0 and
1535 * cycles_wait_est. We will use half of cycles_wait_est.
1537 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1539 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1543 if (fr->nbv->emulateGpu())
1545 // NOTE: emulation kernel is not included in the balancing region,
1546 // but emulation mode does not target performance anyway
1547 wallcycle_start_nocount(wcycle, ewcFORCE);
1548 do_nb_verlet(fr, ic, enerd, forceFlags, Nbnxm::InteractionLocality::Local,
1549 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1550 step, nrnb, wcycle);
1551 wallcycle_stop(wcycle, ewcFORCE);
1554 /* Do the nonbonded GPU (or emulation) force buffer reduction
1555 * on the non-alternating path. */
1556 if (bUseOrEmulGPU && !alternateGpuWait)
1558 gmx::ArrayRef<gmx::RVec> force = forceOut.forceWithShiftForces().force();
1559 rvec *f = as_rvec_array(force.data());
1561 // TODO: move these steps as early as possible:
1562 // - CPU f H2D should be as soon as all CPU-side forces are done
1563 // - wait for force reduction does not need to block host (at least not here, it's sufficient to wait
1564 // before the next CPU task that consumes the forces: vsite spread or update)
1566 if (useGpuFBufOps == BufferOpsUseGpu::True && (haveCpuForces || DOMAINDECOMP(cr)))
1568 nbv->launch_copy_f_to_gpu(f, Nbnxm::AtomLocality::Local);
1570 // flag to specify if forces should be accumulated in force
1571 // buffer ops. For local part, this depends on whether CPU
1572 // forces are present, or if DD is active (in which case the
1573 // halo exchange has resulted in contributions from the
1575 bool accumulateForce = (useGpuFBufOps == BufferOpsUseGpu::True) &&
1576 (haveCpuForces || DOMAINDECOMP(cr));
1577 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
1578 force, pme_gpu_get_device_f(fr->pmedata),
1579 pme_gpu_get_f_ready_synchronizer(fr->pmedata),
1580 useGpuFBufOps, useGpuFPmeReduction, accumulateForce);
1581 if (useGpuFBufOps == BufferOpsUseGpu::True)
1583 nbv->launch_copy_f_from_gpu(f, Nbnxm::AtomLocality::Local);
1584 nbv->wait_for_gpu_force_reduction(Nbnxm::AtomLocality::Local);
1588 launchGpuEndOfStepTasks(nbv, fr->gpuBonded, fr->pmedata, enerd,
1594 if (DOMAINDECOMP(cr))
1596 dd_force_flop_stop(cr->dd, nrnb);
1599 if (forceFlags.computeForces)
1601 rvec *f = as_rvec_array(forceOut.forceWithShiftForces().force().data());
1603 /* If we have NoVirSum forces, but we do not calculate the virial,
1604 * we sum fr->f_novirsum=forceOut.f later.
1606 if (vsite && !(fr->haveDirectVirialContributions && !forceFlags.computeVirial))
1608 rvec *fshift = as_rvec_array(forceOut.forceWithShiftForces().shiftForces().data());
1609 spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), f, fshift, FALSE, nullptr, nrnb,
1610 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1613 if (forceFlags.computeVirial)
1615 /* Calculation of the virial must be done after vsites! */
1616 calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()),
1617 forceOut.forceWithShiftForces(),
1618 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1622 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1624 /* In case of node-splitting, the PP nodes receive the long-range
1625 * forces, virial and energy from the PME nodes here.
1627 pme_receive_force_ener(cr, &forceOut.forceWithVirial(), enerd, wcycle);
1630 if (forceFlags.computeForces)
1632 post_process_forces(cr, step, nrnb, wcycle,
1633 top, box, as_rvec_array(x.unpaddedArrayRef().data()), &forceOut,
1634 vir_force, mdatoms, graph, fr, vsite,
1638 if (forceFlags.computeEnergy)
1640 /* Sum the potential energy terms from group contributions */
1641 sum_epot(&(enerd->grpp), enerd->term);
1643 if (!EI_TPI(inputrec->eI))
1645 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1649 /* In case we don't have constraints and are using GPUs, the next balancing
1650 * region starts here.
1651 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1652 * virial calculation and COM pulling, is not thus not included in
1653 * the balance timing, which is ok as most tasks do communication.
1655 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);