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.
50 #include "gromacs/awh/awh.h"
51 #include "gromacs/domdec/dlbtiming.h"
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/domdec/partition.h"
55 #include "gromacs/essentialdynamics/edsam.h"
56 #include "gromacs/ewald/pme.h"
57 #include "gromacs/gmxlib/chargegroup.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/gmxlib/nrnb.h"
60 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
61 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
62 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
63 #include "gromacs/gpu_utils/gpu_utils.h"
64 #include "gromacs/imd/imd.h"
65 #include "gromacs/listed_forces/bonded.h"
66 #include "gromacs/listed_forces/disre.h"
67 #include "gromacs/listed_forces/gpubonded.h"
68 #include "gromacs/listed_forces/listed_forces.h"
69 #include "gromacs/listed_forces/manage_threading.h"
70 #include "gromacs/listed_forces/orires.h"
71 #include "gromacs/math/arrayrefwithpadding.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/units.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/math/vecdump.h"
76 #include "gromacs/mdlib/calcmu.h"
77 #include "gromacs/mdlib/calcvir.h"
78 #include "gromacs/mdlib/constr.h"
79 #include "gromacs/mdlib/force.h"
80 #include "gromacs/mdlib/forcerec.h"
81 #include "gromacs/mdlib/gmx_omp_nthreads.h"
82 #include "gromacs/mdlib/ppforceworkload.h"
83 #include "gromacs/mdlib/qmmm.h"
84 #include "gromacs/mdlib/update.h"
85 #include "gromacs/mdtypes/commrec.h"
86 #include "gromacs/mdtypes/enerdata.h"
87 #include "gromacs/mdtypes/forceoutput.h"
88 #include "gromacs/mdtypes/iforceprovider.h"
89 #include "gromacs/mdtypes/inputrec.h"
90 #include "gromacs/mdtypes/md_enums.h"
91 #include "gromacs/mdtypes/state.h"
92 #include "gromacs/nbnxm/atomdata.h"
93 #include "gromacs/nbnxm/gpu_data_mgmt.h"
94 #include "gromacs/nbnxm/nbnxm.h"
95 #include "gromacs/pbcutil/ishift.h"
96 #include "gromacs/pbcutil/mshift.h"
97 #include "gromacs/pbcutil/pbc.h"
98 #include "gromacs/pulling/pull.h"
99 #include "gromacs/pulling/pull_rotation.h"
100 #include "gromacs/timing/cyclecounter.h"
101 #include "gromacs/timing/gpu_timing.h"
102 #include "gromacs/timing/wallcycle.h"
103 #include "gromacs/timing/wallcyclereporting.h"
104 #include "gromacs/timing/walltime_accounting.h"
105 #include "gromacs/topology/topology.h"
106 #include "gromacs/utility/arrayref.h"
107 #include "gromacs/utility/basedefinitions.h"
108 #include "gromacs/utility/cstringutil.h"
109 #include "gromacs/utility/exceptions.h"
110 #include "gromacs/utility/fatalerror.h"
111 #include "gromacs/utility/gmxassert.h"
112 #include "gromacs/utility/gmxmpi.h"
113 #include "gromacs/utility/logger.h"
114 #include "gromacs/utility/smalloc.h"
115 #include "gromacs/utility/strconvert.h"
116 #include "gromacs/utility/sysinfo.h"
118 // TODO: this environment variable allows us to verify before release
119 // that on less common architectures the total cost of polling is not larger than
120 // a blocking wait (so polling does not introduce overhead when the static
121 // PME-first ordering would suffice).
122 static const bool c_disableAlternatingWait = (getenv("GMX_DISABLE_ALTERNATING_GPU_WAIT") != nullptr);
125 static void sum_forces(rvec f[], gmx::ArrayRef<const gmx::RVec> forceToAdd)
127 const int end = forceToAdd.size();
129 int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
130 #pragma omp parallel for num_threads(nt) schedule(static)
131 for (int i = 0; i < end; i++)
133 rvec_inc(f[i], forceToAdd[i]);
137 static void calc_virial(int start, int homenr, const rvec x[], const rvec f[],
138 tensor vir_part, const t_graph *graph, const matrix box,
139 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
141 /* The short-range virial from surrounding boxes */
142 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
143 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
145 /* Calculate partial virial, for local atoms only, based on short range.
146 * Total virial is computed in global_stat, called from do_md
148 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
149 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
153 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
157 static void pull_potential_wrapper(const t_commrec *cr,
158 const t_inputrec *ir,
159 const matrix box, gmx::ArrayRef<const gmx::RVec> x,
160 gmx::ForceWithVirial *force,
161 const t_mdatoms *mdatoms,
162 gmx_enerdata_t *enerd,
165 gmx_wallcycle_t wcycle)
170 /* Calculate the center of mass forces, this requires communication,
171 * which is why pull_potential is called close to other communication.
173 wallcycle_start(wcycle, ewcPULLPOT);
174 set_pbc(&pbc, ir->ePBC, box);
176 enerd->term[F_COM_PULL] +=
177 pull_potential(ir->pull_work, mdatoms, &pbc,
178 cr, t, lambda[efptRESTRAINT], as_rvec_array(x.data()), force, &dvdl);
179 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
180 wallcycle_stop(wcycle, ewcPULLPOT);
183 static void pme_receive_force_ener(const t_commrec *cr,
184 gmx::ForceWithVirial *forceWithVirial,
185 gmx_enerdata_t *enerd,
186 gmx_wallcycle_t wcycle)
188 real e_q, e_lj, dvdl_q, dvdl_lj;
189 float cycles_ppdpme, cycles_seppme;
191 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
192 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
194 /* In case of node-splitting, the PP nodes receive the long-range
195 * forces, virial and energy from the PME nodes here.
197 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
200 gmx_pme_receive_f(cr, forceWithVirial, &e_q, &e_lj, &dvdl_q, &dvdl_lj,
202 enerd->term[F_COUL_RECIP] += e_q;
203 enerd->term[F_LJ_RECIP] += e_lj;
204 enerd->dvdl_lin[efptCOUL] += dvdl_q;
205 enerd->dvdl_lin[efptVDW] += dvdl_lj;
209 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
211 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
214 static void print_large_forces(FILE *fp,
222 real force2Tolerance = gmx::square(forceTolerance);
223 gmx::index numNonFinite = 0;
224 for (int i = 0; i < md->homenr; i++)
226 real force2 = norm2(f[i]);
227 bool nonFinite = !std::isfinite(force2);
228 if (force2 >= force2Tolerance || nonFinite)
230 fprintf(fp, "step %" PRId64 " atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
232 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
239 if (numNonFinite > 0)
241 /* Note that with MPI this fatal call on one rank might interrupt
242 * the printing on other ranks. But we can only avoid that with
243 * an expensive MPI barrier that we would need at each step.
245 gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
249 static void post_process_forces(const t_commrec *cr,
252 gmx_wallcycle_t wcycle,
253 const gmx_localtop_t *top,
257 gmx::ForceWithVirial *forceWithVirial,
259 const t_mdatoms *mdatoms,
260 const t_graph *graph,
261 const t_forcerec *fr,
262 const gmx_vsite_t *vsite,
265 if (fr->haveDirectVirialContributions)
267 rvec *fDirectVir = as_rvec_array(forceWithVirial->force_.data());
271 /* Spread the mesh force on virtual sites to the other particles...
272 * This is parallellized. MPI communication is performed
273 * if the constructing atoms aren't local.
275 matrix virial = { { 0 } };
276 spread_vsite_f(vsite, x, fDirectVir, nullptr,
277 (flags & GMX_FORCE_VIRIAL) != 0, virial,
279 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
280 forceWithVirial->addVirialContribution(virial);
283 if (flags & GMX_FORCE_VIRIAL)
285 /* Now add the forces, this is local */
286 sum_forces(f, forceWithVirial->force_);
288 /* Add the direct virial contributions */
289 GMX_ASSERT(forceWithVirial->computeVirial_, "forceWithVirial should request virial computation when we request the virial");
290 m_add(vir_force, forceWithVirial->getVirial(), vir_force);
294 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
299 if (fr->print_force >= 0)
301 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
305 static void do_nb_verlet(t_forcerec *fr,
306 const interaction_const_t *ic,
307 gmx_enerdata_t *enerd,
309 const Nbnxm::InteractionLocality ilocality,
313 gmx_wallcycle_t wcycle)
315 if (!(flags & GMX_FORCE_NONBONDED))
317 /* skip non-bonded calculation */
321 nonbonded_verlet_t *nbv = fr->nbv.get();
323 /* GPU kernel launch overhead is already timed separately */
324 if (fr->cutoff_scheme != ecutsVERLET)
326 gmx_incons("Invalid cut-off scheme passed!");
331 /* When dynamic pair-list pruning is requested, we need to prune
332 * at nstlistPrune steps.
334 if (nbv->pairlistSets().isDynamicPruningStepCpu(step))
336 /* Prune the pair-list beyond fr->ic->rlistPrune using
337 * the current coordinates of the atoms.
339 wallcycle_sub_start(wcycle, ewcsNONBONDED_PRUNING);
340 nbv->dispatchPruneKernelCpu(ilocality, fr->shift_vec);
341 wallcycle_sub_stop(wcycle, ewcsNONBONDED_PRUNING);
344 wallcycle_sub_start(wcycle, ewcsNONBONDED);
347 nbv->dispatchNonbondedKernel(ilocality, *ic, flags, clearF, fr, enerd, nrnb);
351 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
355 gmx_bool use_GPU(const nonbonded_verlet_t *nbv)
357 return nbv != nullptr && nbv->useGpu();
360 static inline void clear_rvecs_omp(int n, rvec v[])
362 int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, n);
364 /* Note that we would like to avoid this conditional by putting it
365 * into the omp pragma instead, but then we still take the full
366 * omp parallel for overhead (at least with gcc5).
370 for (int i = 0; i < n; i++)
377 #pragma omp parallel for num_threads(nth) schedule(static)
378 for (int i = 0; i < n; i++)
385 /*! \brief Return an estimate of the average kinetic energy or 0 when unreliable
387 * \param groupOptions Group options, containing T-coupling options
389 static real averageKineticEnergyEstimate(const t_grpopts &groupOptions)
391 real nrdfCoupled = 0;
392 real nrdfUncoupled = 0;
393 real kineticEnergy = 0;
394 for (int g = 0; g < groupOptions.ngtc; g++)
396 if (groupOptions.tau_t[g] >= 0)
398 nrdfCoupled += groupOptions.nrdf[g];
399 kineticEnergy += groupOptions.nrdf[g]*0.5*groupOptions.ref_t[g]*BOLTZ;
403 nrdfUncoupled += groupOptions.nrdf[g];
407 /* This conditional with > also catches nrdf=0 */
408 if (nrdfCoupled > nrdfUncoupled)
410 return kineticEnergy*(nrdfCoupled + nrdfUncoupled)/nrdfCoupled;
418 /*! \brief This routine checks that the potential energy is finite.
420 * Always checks that the potential energy is finite. If step equals
421 * inputrec.init_step also checks that the magnitude of the potential energy
422 * is reasonable. Terminates with a fatal error when a check fails.
423 * Note that passing this check does not guarantee finite forces,
424 * since those use slightly different arithmetics. But in most cases
425 * there is just a narrow coordinate range where forces are not finite
426 * and energies are finite.
428 * \param[in] step The step number, used for checking and printing
429 * \param[in] enerd The energy data; the non-bonded group energies need to be added to enerd.term[F_EPOT] before calling this routine
430 * \param[in] inputrec The input record
432 static void checkPotentialEnergyValidity(int64_t step,
433 const gmx_enerdata_t &enerd,
434 const t_inputrec &inputrec)
436 /* Threshold valid for comparing absolute potential energy against
437 * the kinetic energy. Normally one should not consider absolute
438 * potential energy values, but with a factor of one million
439 * we should never get false positives.
441 constexpr real c_thresholdFactor = 1e6;
443 bool energyIsNotFinite = !std::isfinite(enerd.term[F_EPOT]);
444 real averageKineticEnergy = 0;
445 /* We only check for large potential energy at the initial step,
446 * because that is by far the most likely step for this too occur
447 * and because computing the average kinetic energy is not free.
448 * Note: nstcalcenergy >> 1 often does not allow to catch large energies
449 * before they become NaN.
451 if (step == inputrec.init_step && EI_DYNAMICS(inputrec.eI))
453 averageKineticEnergy = averageKineticEnergyEstimate(inputrec.opts);
456 if (energyIsNotFinite || (averageKineticEnergy > 0 &&
457 enerd.term[F_EPOT] > c_thresholdFactor*averageKineticEnergy))
459 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.",
462 energyIsNotFinite ? "not finite" : "extremely high",
464 enerd.term[F_COUL_SR],
465 energyIsNotFinite ? "non-finite" : "very high",
466 energyIsNotFinite ? " or Nan" : "");
470 /*! \brief Return true if there are special forces computed this step.
472 * The conditionals exactly correspond to those in computeSpecialForces().
475 haveSpecialForces(const t_inputrec *inputrec,
476 ForceProviders *forceProviders,
480 const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
483 ((computeForces && forceProviders->hasForceProvider()) || // forceProviders
484 (inputrec->bPull && pull_have_potential(inputrec->pull_work)) || // pull
485 inputrec->bRot || // enforced rotation
486 (ed != nullptr) || // flooding
487 (inputrec->bIMD && computeForces)); // IMD
490 /*! \brief Compute forces and/or energies for special algorithms
492 * The intention is to collect all calls to algorithms that compute
493 * forces on local atoms only and that do not contribute to the local
494 * virial sum (but add their virial contribution separately).
495 * Eventually these should likely all become ForceProviders.
496 * Within this function the intention is to have algorithms that do
497 * global communication at the end, so global barriers within the MD loop
498 * are as close together as possible.
500 * \param[in] fplog The log file
501 * \param[in] cr The communication record
502 * \param[in] inputrec The input record
503 * \param[in] awh The Awh module (nullptr if none in use).
504 * \param[in] enforcedRotation Enforced rotation module.
505 * \param[in] step The current MD step
506 * \param[in] t The current time
507 * \param[in,out] wcycle Wallcycle accounting struct
508 * \param[in,out] forceProviders Pointer to a list of force providers
509 * \param[in] box The unit cell
510 * \param[in] x The coordinates
511 * \param[in] mdatoms Per atom properties
512 * \param[in] lambda Array of free-energy lambda values
513 * \param[in] forceFlags Flags that tell whether we should compute forces/energies/virial
514 * \param[in,out] forceWithVirial Force and virial buffers
515 * \param[in,out] enerd Energy buffer
516 * \param[in,out] ed Essential dynamics pointer
517 * \param[in] bNS Tells if we did neighbor searching this step, used for ED sampling
519 * \todo Remove bNS, which is used incorrectly.
520 * \todo Convert all other algorithms called here to ForceProviders.
523 computeSpecialForces(FILE *fplog,
525 const t_inputrec *inputrec,
527 gmx_enfrot *enforcedRotation,
530 gmx_wallcycle_t wcycle,
531 ForceProviders *forceProviders,
533 gmx::ArrayRef<const gmx::RVec> x,
534 const t_mdatoms *mdatoms,
537 gmx::ForceWithVirial *forceWithVirial,
538 gmx_enerdata_t *enerd,
542 const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
544 /* NOTE: Currently all ForceProviders only provide forces.
545 * When they also provide energies, remove this conditional.
549 gmx::ForceProviderInput forceProviderInput(x, *mdatoms, t, box, *cr);
550 gmx::ForceProviderOutput forceProviderOutput(forceWithVirial, enerd);
552 /* Collect forces from modules */
553 forceProviders->calculateForces(forceProviderInput, &forceProviderOutput);
556 if (inputrec->bPull && pull_have_potential(inputrec->pull_work))
558 pull_potential_wrapper(cr, inputrec, box, x,
560 mdatoms, enerd, lambda, t,
565 enerd->term[F_COM_PULL] +=
566 awh->applyBiasForcesAndUpdateBias(inputrec->ePBC, *mdatoms, box,
568 t, step, wcycle, fplog);
572 rvec *f = as_rvec_array(forceWithVirial->force_.data());
574 /* Add the forces from enforced rotation potentials (if any) */
577 wallcycle_start(wcycle, ewcROTadd);
578 enerd->term[F_COM_PULL] += add_rot_forces(enforcedRotation, f, cr, step, t);
579 wallcycle_stop(wcycle, ewcROTadd);
584 /* Note that since init_edsam() is called after the initialization
585 * of forcerec, edsam doesn't request the noVirSum force buffer.
586 * Thus if no other algorithm (e.g. PME) requires it, the forces
587 * here will contribute to the virial.
589 do_flood(cr, inputrec, as_rvec_array(x.data()), f, ed, box, step, bNS);
592 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
593 if (inputrec->bIMD && computeForces)
595 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
599 /*! \brief Launch the prepare_step and spread stages of PME GPU.
601 * \param[in] pmedata The PME structure
602 * \param[in] box The box matrix
603 * \param[in] x Coordinate array
604 * \param[in] flags Force flags
605 * \param[in] pmeFlags PME flags
606 * \param[in] wcycle The wallcycle structure
608 static inline void launchPmeGpuSpread(gmx_pme_t *pmedata,
613 gmx_wallcycle_t wcycle)
615 pme_gpu_prepare_computation(pmedata, (flags & GMX_FORCE_DYNAMICBOX) != 0, box, wcycle, pmeFlags);
616 pme_gpu_launch_spread(pmedata, x, wcycle);
619 /*! \brief Launch the FFT and gather stages of PME GPU
621 * This function only implements setting the output forces (no accumulation).
623 * \param[in] pmedata The PME structure
624 * \param[in] wcycle The wallcycle structure
626 static void launchPmeGpuFftAndGather(gmx_pme_t *pmedata,
627 gmx_wallcycle_t wcycle)
629 pme_gpu_launch_complex_transforms(pmedata, wcycle);
630 pme_gpu_launch_gather(pmedata, wcycle, PmeForceOutputHandling::Set);
634 * Polling wait for either of the PME or nonbonded GPU tasks.
636 * Instead of a static order in waiting for GPU tasks, this function
637 * polls checking which of the two tasks completes first, and does the
638 * associated force buffer reduction overlapped with the other task.
639 * By doing that, unlike static scheduling order, it can always overlap
640 * one of the reductions, regardless of the GPU task completion order.
642 * \param[in] nbv Nonbonded verlet structure
643 * \param[in,out] pmedata PME module data
644 * \param[in,out] force Force array to reduce task outputs into.
645 * \param[in,out] forceWithVirial Force and virial buffers
646 * \param[in,out] fshift Shift force output vector results are reduced into
647 * \param[in,out] enerd Energy data structure results are reduced into
648 * \param[in] flags Force flags
649 * \param[in] pmeFlags PME flags
650 * \param[in] haveOtherWork Tells whether there is other work than non-bonded in the stream(s)
651 * \param[in] wcycle The wallcycle structure
653 static void alternatePmeNbGpuWaitReduce(nonbonded_verlet_t *nbv,
655 gmx::ArrayRefWithPadding<gmx::RVec> *force,
656 gmx::ForceWithVirial *forceWithVirial,
658 gmx_enerdata_t *enerd,
662 gmx_wallcycle_t wcycle)
664 bool isPmeGpuDone = false;
665 bool isNbGpuDone = false;
668 gmx::ArrayRef<const gmx::RVec> pmeGpuForces;
670 while (!isPmeGpuDone || !isNbGpuDone)
674 GpuTaskCompletion completionType = (isNbGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
675 isPmeGpuDone = pme_gpu_try_finish_task(pmedata, pmeFlags, wcycle, forceWithVirial, enerd, completionType);
680 GpuTaskCompletion completionType = (isPmeGpuDone) ? GpuTaskCompletion::Wait : GpuTaskCompletion::Check;
681 wallcycle_start_nocount(wcycle, ewcWAIT_GPU_NB_L);
682 isNbGpuDone = Nbnxm::gpu_try_finish_task(nbv->gpu_nbv,
684 Nbnxm::AtomLocality::Local,
686 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
687 fshift, completionType);
688 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
689 // To get the call count right, when the task finished we
690 // issue a start/stop.
691 // TODO: move the ewcWAIT_GPU_NB_L cycle counting into nbnxn_gpu_try_finish_task()
692 // and ewcNB_XF_BUF_OPS counting into nbnxn_atomdata_add_nbat_f_to_f().
695 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
696 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
698 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
699 as_rvec_array(force->unpaddedArrayRef().data()), wcycle);
705 /*! \brief Hack structure with force ouput buffers for do_force */
709 ForceOutputs(rvec *f, gmx::ForceWithVirial const forceWithVirial) :
711 forceWithVirial(forceWithVirial) {}
713 //! Force output buffer used by legacy modules
715 //! Force with direct virial contribution (if there are any)
716 gmx::ForceWithVirial forceWithVirial;
719 /*! \brief Set up the different force buffers; also does clearing.
721 * \param[in] fr force record pointer
722 * \param[in] inputrec input record
723 * \param[in] force force array
724 * \param[in] bDoForces True if force are computed this step
725 * \param[in] doVirial True if virial is computed this step
726 * \param[out] wcycle wallcycle recording structure
728 * \returns Cleared force output structure
731 setupForceOutputs(const t_forcerec *fr,
732 const t_inputrec &inputrec,
733 gmx::ArrayRefWithPadding<gmx::RVec> force,
734 const bool bDoForces,
736 gmx_wallcycle_t wcycle)
738 wallcycle_sub_start(wcycle, ewcsCLEAR_FORCE_BUFFER);
740 /* Temporary solution until all routines take PaddedRVecVector */
741 rvec *const f = as_rvec_array(force.unpaddedArrayRef().data());
744 /* Clear the short- and long-range forces */
745 clear_rvecs_omp(fr->natoms_force_constr, f);
748 /* If we need to compute the virial, we might need a separate
749 * force buffer for algorithms for which the virial is calculated
750 * directly, such as PME. Otherwise, forceWithVirial uses the
751 * the same force (f in legacy calls) buffer as other algorithms.
753 const bool useSeparateForceWithVirialBuffer = (bDoForces && (doVirial && fr->haveDirectVirialContributions));
756 /* forceWithVirial uses the local atom range only */
757 gmx::ForceWithVirial forceWithVirial (useSeparateForceWithVirialBuffer ?
758 *fr->forceBufferForDirectVirialContributions : force.unpaddedArrayRef(),
761 if (useSeparateForceWithVirialBuffer)
763 /* TODO: update comment
764 * We only compute forces on local atoms. Note that vsites can
765 * spread to non-local atoms, but that part of the buffer is
766 * cleared separately in the vsite spreading code.
768 clear_rvecs_omp(forceWithVirial.force_.size(), as_rvec_array(forceWithVirial.force_.data()));
771 if (inputrec.bPull && pull_have_constraint(inputrec.pull_work))
773 clear_pull_forces(inputrec.pull_work);
776 wallcycle_sub_stop(wcycle, ewcsCLEAR_FORCE_BUFFER);
778 return ForceOutputs(f, forceWithVirial);
782 /*! \brief Set up flags that indicate what type of work is there to compute.
784 * Currently we only update it at search steps,
785 * but some properties may change more frequently (e.g. virial/non-virial step),
786 * so when including those either the frequency of update (per-step) or the scope
787 * of a flag will change (i.e. a set of flags for nstlist steps).
791 setupForceWorkload(gmx::PpForceWorkload *forceWork,
792 const t_inputrec *inputrec,
793 const t_forcerec *fr,
800 forceWork->haveSpecialForces = haveSpecialForces(inputrec, fr->forceProviders, forceFlags, ed);
801 forceWork->haveCpuBondedWork = haveCpuBondeds(*fr);
802 forceWork->haveGpuBondedWork = ((fr->gpuBonded != nullptr) && fr->gpuBonded->haveInteractions());
803 forceWork->haveRestraintsWork = havePositionRestraints(idef, *fcd);
804 forceWork->haveCpuListedForceWork = haveCpuListedForces(*fr, idef, *fcd);
807 static void do_force_cutsVERLET(FILE *fplog,
809 const gmx_multisim_t *ms,
810 const t_inputrec *inputrec,
812 gmx_enfrot *enforcedRotation,
815 gmx_wallcycle_t wcycle,
816 const gmx_localtop_t *top,
817 const gmx_groups_t * /* groups */,
818 matrix box, gmx::ArrayRefWithPadding<gmx::RVec> x,
820 gmx::ArrayRefWithPadding<gmx::RVec> force,
822 const t_mdatoms *mdatoms,
823 gmx_enerdata_t *enerd, t_fcdata *fcd,
827 gmx::PpForceWorkload *ppForceWorkload,
828 interaction_const_t *ic,
829 const gmx_vsite_t *vsite,
834 const DDBalanceRegionHandler &ddBalanceRegionHandler)
838 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
839 gmx_bool bDoForces, bUseGPU, bUseOrEmulGPU;
840 rvec vzero, box_diag;
841 float cycles_pme, cycles_wait_gpu;
842 nonbonded_verlet_t *nbv = fr->nbv.get();
844 bStateChanged = ((flags & GMX_FORCE_STATECHANGED) != 0);
845 bNS = ((flags & GMX_FORCE_NS) != 0);
846 bFillGrid = (bNS && bStateChanged);
847 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
848 bDoForces = ((flags & GMX_FORCE_FORCES) != 0);
849 bUseGPU = fr->nbv->useGpu();
850 bUseOrEmulGPU = bUseGPU || fr->nbv->emulateGpu();
852 const auto pmeRunMode = fr->pmedata ? pme_run_mode(fr->pmedata) : PmeRunMode::CPU;
853 // TODO slim this conditional down - inputrec and duty checks should mean the same in proper code!
854 const bool useGpuPme = EEL_PME(fr->ic->eeltype) && thisRankHasDuty(cr, DUTY_PME) &&
855 ((pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed));
856 const int pmeFlags = GMX_PME_SPREAD | GMX_PME_SOLVE |
857 ((flags & GMX_FORCE_VIRIAL) ? GMX_PME_CALC_ENER_VIR : 0) |
858 ((flags & GMX_FORCE_ENERGY) ? GMX_PME_CALC_ENER_VIR : 0) |
859 ((flags & GMX_FORCE_FORCES) ? GMX_PME_CALC_F : 0);
861 /* At a search step we need to start the first balancing region
862 * somewhere early inside the step after communication during domain
863 * decomposition (and not during the previous step as usual).
867 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::yes);
873 const int homenr = mdatoms->homenr;
875 clear_mat(vir_force);
877 if (DOMAINDECOMP(cr))
879 cg1 = cr->dd->globalAtomGroupIndices.size();
892 update_forcerec(fr, box);
894 if (inputrecNeedMutot(inputrec))
896 /* Calculate total (local) dipole moment in a temporary common array.
897 * This makes it possible to sum them over nodes faster.
899 calc_mu(start, homenr,
900 x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
905 if (fr->ePBC != epbcNONE)
907 /* Compute shift vectors every step,
908 * because of pressure coupling or box deformation!
910 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
912 calc_shifts(box, fr->shift_vec);
917 put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr));
918 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
920 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
922 unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
926 nbnxn_atomdata_copy_shiftvec((flags & GMX_FORCE_DYNAMICBOX) != 0,
927 fr->shift_vec, nbv->nbat.get());
930 if (!thisRankHasDuty(cr, DUTY_PME))
932 /* Send particle coordinates to the pme nodes.
933 * Since this is only implemented for domain decomposition
934 * and domain decomposition does not use the graph,
935 * we do not need to worry about shifting.
937 gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
938 lambda[efptCOUL], lambda[efptVDW],
939 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)) != 0,
946 launchPmeGpuSpread(fr->pmedata, box, as_rvec_array(x.unpaddedArrayRef().data()), flags, pmeFlags, wcycle);
949 /* do gridding for pair search */
952 if (graph && bStateChanged)
954 /* Calculate intramolecular shift vectors to make molecules whole */
955 mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
959 box_diag[XX] = box[XX][XX];
960 box_diag[YY] = box[YY][YY];
961 box_diag[ZZ] = box[ZZ][ZZ];
963 wallcycle_start(wcycle, ewcNS);
964 if (!DOMAINDECOMP(cr))
966 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
967 nbnxn_put_on_grid(nbv, box,
969 nullptr, 0, mdatoms->homenr, -1,
970 fr->cginfo, x.unpaddedArrayRef(),
972 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
976 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
977 nbnxn_put_on_grid_nonlocal(nbv, domdec_zones(cr->dd),
978 fr->cginfo, x.unpaddedArrayRef());
979 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
982 nbv->setAtomProperties(*mdatoms, *fr->cginfo);
984 wallcycle_stop(wcycle, ewcNS);
986 /* initialize the GPU nbnxm atom data and bonded data structures */
989 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
991 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
992 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
993 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
997 /* Now we put all atoms on the grid, we can assign bonded
998 * interactions to the GPU, where the grid order is
999 * needed. Also the xq, f and fshift device buffers have
1000 * been reallocated if needed, so the bonded code can
1001 * learn about them. */
1002 // TODO the xq, f, and fshift buffers are now shared
1003 // resources, so they should be maintained by a
1004 // higher-level object than the nb module.
1005 fr->gpuBonded->updateInteractionListsAndDeviceBuffers(nbv->getGridIndices(),
1007 Nbnxm::gpu_get_xq(nbv->gpu_nbv),
1008 Nbnxm::gpu_get_f(nbv->gpu_nbv),
1009 Nbnxm::gpu_get_fshift(nbv->gpu_nbv));
1011 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1014 // Need to run after the GPU-offload bonded interaction lists
1015 // are set up to be able to determine whether there is bonded work.
1016 setupForceWorkload(ppForceWorkload,
1025 /* do local pair search */
1028 // TODO: fuse this branch with the above bNS block
1029 wallcycle_start_nocount(wcycle, ewcNS);
1030 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1031 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1032 nbv->constructPairlist(Nbnxm::InteractionLocality::Local,
1033 &top->excls, step, nrnb);
1034 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1035 wallcycle_stop(wcycle, ewcNS);
1039 nbv->setCoordinates(Nbnxm::AtomLocality::Local, false,
1040 x.unpaddedArrayRef(), wcycle);
1045 ddBalanceRegionHandler.openBeforeForceComputationGpu();
1047 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1049 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1050 Nbnxm::gpu_upload_shiftvec(nbv->gpu_nbv, nbv->nbat.get());
1051 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1052 Nbnxm::AtomLocality::Local,
1053 ppForceWorkload->haveGpuBondedWork);
1054 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1056 // bonded work not split into separate local and non-local, so with DD
1057 // we can only launch the kernel after non-local coordinates have been received.
1058 if (ppForceWorkload->haveGpuBondedWork && !havePPDomainDecomposition(cr))
1060 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1061 fr->gpuBonded->launchKernels(fr, flags, box);
1062 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1065 /* launch local nonbonded work on GPU */
1066 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1067 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFNo,
1068 step, nrnb, wcycle);
1069 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1070 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1075 // In PME GPU and mixed mode we launch FFT / gather after the
1076 // X copy/transform to allow overlap as well as after the GPU NB
1077 // launch to avoid FFT launch overhead hijacking the CPU and delaying
1078 // the nonbonded kernel.
1079 launchPmeGpuFftAndGather(fr->pmedata, wcycle);
1082 /* Communicate coordinates and sum dipole if necessary +
1083 do non-local pair search */
1084 if (havePPDomainDecomposition(cr))
1088 // TODO: fuse this branch with the above large bNS block
1089 wallcycle_start_nocount(wcycle, ewcNS);
1090 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1091 /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
1092 nbv->constructPairlist(Nbnxm::InteractionLocality::NonLocal,
1093 &top->excls, step, nrnb);
1094 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1095 wallcycle_stop(wcycle, ewcNS);
1099 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1101 nbv->setCoordinates(Nbnxm::AtomLocality::NonLocal, false,
1102 x.unpaddedArrayRef(), wcycle);
1107 wallcycle_start(wcycle, ewcLAUNCH_GPU);
1109 /* launch non-local nonbonded tasks on GPU */
1110 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1111 Nbnxm::gpu_copy_xq_to_gpu(nbv->gpu_nbv, nbv->nbat.get(),
1112 Nbnxm::AtomLocality::NonLocal,
1113 ppForceWorkload->haveGpuBondedWork);
1114 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1116 if (ppForceWorkload->haveGpuBondedWork)
1118 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_BONDED);
1119 fr->gpuBonded->launchKernels(fr, flags, box);
1120 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1123 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1124 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1125 step, nrnb, wcycle);
1126 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1128 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1134 /* launch D2H copy-back F */
1135 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1136 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1137 if (havePPDomainDecomposition(cr))
1139 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1140 flags, Nbnxm::AtomLocality::NonLocal, ppForceWorkload->haveGpuBondedWork);
1142 Nbnxm::gpu_launch_cpyback(nbv->gpu_nbv, nbv->nbat.get(),
1143 flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork);
1144 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1146 if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1148 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1149 fr->gpuBonded->launchEnergyTransfer();
1150 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1152 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1155 if (bStateChanged && inputrecNeedMutot(inputrec))
1159 gmx_sumd(2*DIM, mu, cr);
1161 ddBalanceRegionHandler.reopenRegionCpu();
1164 for (i = 0; i < 2; i++)
1166 for (j = 0; j < DIM; j++)
1168 fr->mu_tot[i][j] = mu[i*DIM + j];
1172 if (fr->efep == efepNO)
1174 copy_rvec(fr->mu_tot[0], mu_tot);
1178 for (j = 0; j < DIM; j++)
1181 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1182 lambda[efptCOUL]*fr->mu_tot[1][j];
1186 /* Reset energies */
1187 reset_enerdata(enerd);
1188 clear_rvecs(SHIFTS, fr->fshift);
1190 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1192 wallcycle_start(wcycle, ewcPPDURINGPME);
1193 dd_force_flop_start(cr->dd, nrnb);
1198 wallcycle_start(wcycle, ewcROT);
1199 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, bNS);
1200 wallcycle_stop(wcycle, ewcROT);
1203 /* Start the force cycle counter.
1204 * Note that a different counter is used for dynamic load balancing.
1206 wallcycle_start(wcycle, ewcFORCE);
1208 // set up and clear force outputs
1209 struct ForceOutputs forceOut = setupForceOutputs(fr, *inputrec, force, bDoForces, ((flags & GMX_FORCE_VIRIAL) != 0), wcycle);
1211 /* We calculate the non-bonded forces, when done on the CPU, here.
1212 * We do this before calling do_force_lowlevel, because in that
1213 * function, the listed forces are calculated before PME, which
1214 * does communication. With this order, non-bonded and listed
1215 * force calculation imbalance can be balanced out by the domain
1216 * decomposition load balancing.
1221 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local, enbvClearFYes,
1222 step, nrnb, wcycle);
1225 if (fr->efep != efepNO)
1227 /* Calculate the local and non-local free energy interactions here.
1228 * Happens here on the CPU both with and without GPU.
1230 wallcycle_sub_start(wcycle, ewcsNONBONDED);
1231 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::Local,
1232 fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1233 inputrec->fepvals, lambda,
1234 enerd, flags, nrnb);
1236 if (havePPDomainDecomposition(cr))
1238 nbv->dispatchFreeEnergyKernel(Nbnxm::InteractionLocality::NonLocal,
1239 fr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, *mdatoms,
1240 inputrec->fepvals, lambda,
1241 enerd, flags, nrnb);
1243 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
1248 if (havePPDomainDecomposition(cr))
1250 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFNo,
1251 step, nrnb, wcycle);
1254 /* Add all the non-bonded force to the normal force array.
1255 * This can be split into a local and a non-local part when overlapping
1256 * communication with calculation with domain decomposition.
1258 wallcycle_stop(wcycle, ewcFORCE);
1260 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::All, forceOut.f, wcycle);
1262 wallcycle_start_nocount(wcycle, ewcFORCE);
1264 /* If there are multiple fshift output buffers we need to reduce them */
1265 if (flags & GMX_FORCE_VIRIAL)
1267 /* This is not in a subcounter because it takes a
1268 negligible and constant-sized amount of time */
1269 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->nbat.get(),
1274 /* update QMMMrec, if necessary */
1277 update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1280 /* Compute the bonded and non-bonded energies and optionally forces */
1281 do_force_lowlevel(fr, inputrec, &(top->idef),
1282 cr, ms, nrnb, wcycle, mdatoms,
1283 as_rvec_array(x.unpaddedArrayRef().data()), hist, forceOut.f, &forceOut.forceWithVirial, enerd, fcd,
1284 box, inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1286 &cycles_pme, ddBalanceRegionHandler);
1288 wallcycle_stop(wcycle, ewcFORCE);
1290 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1292 fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda,
1293 flags, &forceOut.forceWithVirial, enerd,
1298 /* wait for non-local forces (or calculate in emulation mode) */
1299 if (havePPDomainDecomposition(cr))
1303 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1304 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1305 flags, Nbnxm::AtomLocality::NonLocal,
1306 ppForceWorkload->haveGpuBondedWork,
1307 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1309 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1313 wallcycle_start_nocount(wcycle, ewcFORCE);
1314 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::NonLocal, enbvClearFYes,
1315 step, nrnb, wcycle);
1316 wallcycle_stop(wcycle, ewcFORCE);
1319 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::NonLocal,
1320 forceOut.f, wcycle);
1324 if (havePPDomainDecomposition(cr))
1326 /* We are done with the CPU compute.
1327 * We will now communicate the non-local forces.
1328 * If we use a GPU this will overlap with GPU work, so in that case
1329 * we do not close the DD force balancing region here.
1331 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1335 dd_move_f(cr->dd, force.unpaddedArrayRef(), fr->fshift, wcycle);
1339 // With both nonbonded and PME offloaded a GPU on the same rank, we use
1340 // an alternating wait/reduction scheme.
1341 bool alternateGpuWait = (!c_disableAlternatingWait && useGpuPme && bUseGPU && !DOMAINDECOMP(cr));
1342 if (alternateGpuWait)
1344 alternatePmeNbGpuWaitReduce(fr->nbv.get(), fr->pmedata, &force, &forceOut.forceWithVirial, fr->fshift, enerd,
1345 flags, pmeFlags, ppForceWorkload->haveGpuBondedWork, wcycle);
1348 if (!alternateGpuWait && useGpuPme)
1350 pme_gpu_wait_and_reduce(fr->pmedata, pmeFlags, wcycle, &forceOut.forceWithVirial, enerd);
1353 /* Wait for local GPU NB outputs on the non-alternating wait path */
1354 if (!alternateGpuWait && bUseGPU)
1356 /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1357 * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1358 * but even with a step of 0.1 ms the difference is less than 1%
1361 const float gpuWaitApiOverheadMargin = 2e6f; /* cycles */
1363 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1364 Nbnxm::gpu_wait_finish_task(nbv->gpu_nbv,
1365 flags, Nbnxm::AtomLocality::Local, ppForceWorkload->haveGpuBondedWork,
1366 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1368 float cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1370 if (ddBalanceRegionHandler.useBalancingRegion())
1372 DdBalanceRegionWaitedForGpu waitedForGpu = DdBalanceRegionWaitedForGpu::yes;
1373 if (bDoForces && cycles_tmp <= gpuWaitApiOverheadMargin)
1375 /* We measured few cycles, it could be that the kernel
1376 * and transfer finished earlier and there was no actual
1377 * wait time, only API call overhead.
1378 * Then the actual time could be anywhere between 0 and
1379 * cycles_wait_est. We will use half of cycles_wait_est.
1381 waitedForGpu = DdBalanceRegionWaitedForGpu::no;
1383 ddBalanceRegionHandler.closeAfterForceComputationGpu(cycles_wait_gpu, waitedForGpu);
1387 if (fr->nbv->emulateGpu())
1389 // NOTE: emulation kernel is not included in the balancing region,
1390 // but emulation mode does not target performance anyway
1391 wallcycle_start_nocount(wcycle, ewcFORCE);
1392 do_nb_verlet(fr, ic, enerd, flags, Nbnxm::InteractionLocality::Local,
1393 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1394 step, nrnb, wcycle);
1395 wallcycle_stop(wcycle, ewcFORCE);
1400 pme_gpu_reinit_computation(fr->pmedata, wcycle);
1405 /* now clear the GPU outputs while we finish the step on the CPU */
1406 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1407 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1408 Nbnxm::gpu_clear_outputs(nbv->gpu_nbv, flags);
1410 if (nbv->pairlistSets().isDynamicPruningStepGpu(step))
1412 nbv->dispatchPruneKernelGpu(step);
1414 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_NONBONDED);
1415 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1418 if (ppForceWorkload->haveGpuBondedWork && (flags & GMX_FORCE_ENERGY))
1420 wallcycle_start(wcycle, ewcWAIT_GPU_BONDED);
1421 // in principle this should be included in the DD balancing region,
1422 // but generally it is infrequent so we'll omit it for the sake of
1424 fr->gpuBonded->accumulateEnergyTerms(enerd);
1425 wallcycle_stop(wcycle, ewcWAIT_GPU_BONDED);
1427 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
1428 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_BONDED);
1429 fr->gpuBonded->clearEnergies();
1430 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_BONDED);
1431 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
1434 /* Do the nonbonded GPU (or emulation) force buffer reduction
1435 * on the non-alternating path. */
1436 if (bUseOrEmulGPU && !alternateGpuWait)
1438 nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::Local,
1439 forceOut.f, wcycle);
1441 if (DOMAINDECOMP(cr))
1443 dd_force_flop_stop(cr->dd, nrnb);
1448 /* If we have NoVirSum forces, but we do not calculate the virial,
1449 * we sum fr->f_novirsum=forceOut.f later.
1451 if (vsite && !(fr->haveDirectVirialContributions && !(flags & GMX_FORCE_VIRIAL)))
1453 spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, fr->fshift, FALSE, nullptr, nrnb,
1454 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1457 if (flags & GMX_FORCE_VIRIAL)
1459 /* Calculation of the virial must be done after vsites! */
1460 calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f,
1461 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1465 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1467 /* In case of node-splitting, the PP nodes receive the long-range
1468 * forces, virial and energy from the PME nodes here.
1470 pme_receive_force_ener(cr, &forceOut.forceWithVirial, enerd, wcycle);
1475 post_process_forces(cr, step, nrnb, wcycle,
1476 top, box, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, &forceOut.forceWithVirial,
1477 vir_force, mdatoms, graph, fr, vsite,
1481 if (flags & GMX_FORCE_ENERGY)
1483 /* Sum the potential energy terms from group contributions */
1484 sum_epot(&(enerd->grpp), enerd->term);
1486 if (!EI_TPI(inputrec->eI))
1488 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1493 static void do_force_cutsGROUP(FILE *fplog,
1494 const t_commrec *cr,
1495 const gmx_multisim_t *ms,
1496 const t_inputrec *inputrec,
1498 gmx_enfrot *enforcedRotation,
1501 gmx_wallcycle_t wcycle,
1502 gmx_localtop_t *top,
1503 const gmx_groups_t *groups,
1504 matrix box, gmx::ArrayRefWithPadding<gmx::RVec> x,
1506 gmx::ArrayRefWithPadding<gmx::RVec> force,
1508 const t_mdatoms *mdatoms,
1509 gmx_enerdata_t *enerd,
1514 const gmx_vsite_t *vsite,
1519 const DDBalanceRegionHandler &ddBalanceRegionHandler)
1523 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
1527 const int start = 0;
1528 const int homenr = mdatoms->homenr;
1530 clear_mat(vir_force);
1533 if (DOMAINDECOMP(cr))
1535 cg1 = cr->dd->globalAtomGroupIndices.size();
1546 bStateChanged = ((flags & GMX_FORCE_STATECHANGED) != 0);
1547 bNS = ((flags & GMX_FORCE_NS) != 0);
1548 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1549 bFillGrid = (bNS && bStateChanged);
1550 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1551 bDoForces = ((flags & GMX_FORCE_FORCES) != 0);
1555 update_forcerec(fr, box);
1557 if (inputrecNeedMutot(inputrec))
1559 /* Calculate total (local) dipole moment in a temporary common array.
1560 * This makes it possible to sum them over nodes faster.
1562 calc_mu(start, homenr,
1563 x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1568 if (fr->ePBC != epbcNONE)
1570 /* Compute shift vectors every step,
1571 * because of pressure coupling or box deformation!
1573 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1575 calc_shifts(box, fr->shift_vec);
1580 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1581 &(top->cgs), as_rvec_array(x.unpaddedArrayRef().data()), fr->cg_cm);
1582 inc_nrnb(nrnb, eNR_CGCM, homenr);
1583 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1585 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1587 unshift_self(graph, box, as_rvec_array(x.unpaddedArrayRef().data()));
1592 calc_cgcm(fplog, cg0, cg1, &(top->cgs), as_rvec_array(x.unpaddedArrayRef().data()), fr->cg_cm);
1593 inc_nrnb(nrnb, eNR_CGCM, homenr);
1596 if (bCalcCGCM && gmx_debug_at)
1598 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1602 if (!thisRankHasDuty(cr, DUTY_PME))
1604 /* Send particle coordinates to the pme nodes.
1605 * Since this is only implemented for domain decomposition
1606 * and domain decomposition does not use the graph,
1607 * we do not need to worry about shifting.
1609 gmx_pme_send_coordinates(cr, box, as_rvec_array(x.unpaddedArrayRef().data()),
1610 lambda[efptCOUL], lambda[efptVDW],
1611 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)) != 0,
1614 #endif /* GMX_MPI */
1616 /* Communicate coordinates and sum dipole if necessary */
1617 if (DOMAINDECOMP(cr))
1619 dd_move_x(cr->dd, box, x.unpaddedArrayRef(), wcycle);
1621 /* No GPU support, no move_x overlap, so reopen the balance region here */
1622 ddBalanceRegionHandler.reopenRegionCpu();
1625 if (inputrecNeedMutot(inputrec))
1631 gmx_sumd(2*DIM, mu, cr);
1633 ddBalanceRegionHandler.reopenRegionCpu();
1635 for (i = 0; i < 2; i++)
1637 for (j = 0; j < DIM; j++)
1639 fr->mu_tot[i][j] = mu[i*DIM + j];
1643 if (fr->efep == efepNO)
1645 copy_rvec(fr->mu_tot[0], mu_tot);
1649 for (j = 0; j < DIM; j++)
1652 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1657 /* Reset energies */
1658 reset_enerdata(enerd);
1659 clear_rvecs(SHIFTS, fr->fshift);
1663 wallcycle_start(wcycle, ewcNS);
1665 if (graph && bStateChanged)
1667 /* Calculate intramolecular shift vectors to make molecules whole */
1668 mk_mshift(fplog, graph, fr->ePBC, box, as_rvec_array(x.unpaddedArrayRef().data()));
1671 /* Do the actual neighbour searching */
1673 groups, top, mdatoms,
1674 cr, nrnb, bFillGrid);
1676 wallcycle_stop(wcycle, ewcNS);
1679 if (DOMAINDECOMP(cr) && !thisRankHasDuty(cr, DUTY_PME))
1681 wallcycle_start(wcycle, ewcPPDURINGPME);
1682 dd_force_flop_start(cr->dd, nrnb);
1687 wallcycle_start(wcycle, ewcROT);
1688 do_rotation(cr, enforcedRotation, box, as_rvec_array(x.unpaddedArrayRef().data()), t, step, bNS);
1689 wallcycle_stop(wcycle, ewcROT);
1693 /* Start the force cycle counter.
1694 * Note that a different counter is used for dynamic load balancing.
1696 wallcycle_start(wcycle, ewcFORCE);
1698 // set up and clear force outputs
1699 struct ForceOutputs forceOut = setupForceOutputs(fr, *inputrec, force, bDoForces, ((flags & GMX_FORCE_VIRIAL) != 0), wcycle);
1701 if (inputrec->bPull && pull_have_constraint(inputrec->pull_work))
1703 clear_pull_forces(inputrec->pull_work);
1706 /* update QMMMrec, if necessary */
1709 update_QMMMrec(cr, fr, as_rvec_array(x.unpaddedArrayRef().data()), mdatoms, box);
1712 /* Compute the bonded and non-bonded energies and optionally forces */
1713 do_force_lowlevel(fr, inputrec, &(top->idef),
1714 cr, ms, nrnb, wcycle, mdatoms,
1715 as_rvec_array(x.unpaddedArrayRef().data()), hist, forceOut.f, &forceOut.forceWithVirial, enerd, fcd,
1716 box, inputrec->fepvals, lambda,
1717 graph, &(top->excls), fr->mu_tot,
1719 &cycles_pme, ddBalanceRegionHandler);
1721 wallcycle_stop(wcycle, ewcFORCE);
1723 if (DOMAINDECOMP(cr))
1725 dd_force_flop_stop(cr->dd, nrnb);
1727 ddBalanceRegionHandler.closeAfterForceComputationCpu();
1730 computeSpecialForces(fplog, cr, inputrec, awh, enforcedRotation,
1732 fr->forceProviders, box, x.unpaddedArrayRef(), mdatoms, lambda,
1733 flags, &forceOut.forceWithVirial, enerd,
1738 /* Communicate the forces */
1739 if (DOMAINDECOMP(cr))
1741 dd_move_f(cr->dd, force.unpaddedArrayRef(), fr->fshift, wcycle);
1742 /* Do we need to communicate the separate force array
1743 * for terms that do not contribute to the single sum virial?
1744 * Position restraints and electric fields do not introduce
1745 * inter-cg forces, only full electrostatics methods do.
1746 * When we do not calculate the virial, fr->f_novirsum = forceOut.f,
1747 * so we have already communicated these forces.
1749 if (EEL_FULL(fr->ic->eeltype) && cr->dd->n_intercg_excl &&
1750 (flags & GMX_FORCE_VIRIAL))
1752 dd_move_f(cr->dd, forceOut.forceWithVirial.force_, nullptr, wcycle);
1756 /* If we have NoVirSum forces, but we do not calculate the virial,
1757 * we sum fr->f_novirsum=forceOut.f later.
1759 if (vsite && !(fr->haveDirectVirialContributions && !(flags & GMX_FORCE_VIRIAL)))
1761 spread_vsite_f(vsite, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, fr->fshift, FALSE, nullptr, nrnb,
1762 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr, wcycle);
1765 if (flags & GMX_FORCE_VIRIAL)
1767 /* Calculation of the virial must be done after vsites! */
1768 calc_virial(0, mdatoms->homenr, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f,
1769 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1773 if (PAR(cr) && !thisRankHasDuty(cr, DUTY_PME))
1775 /* In case of node-splitting, the PP nodes receive the long-range
1776 * forces, virial and energy from the PME nodes here.
1778 pme_receive_force_ener(cr, &forceOut.forceWithVirial, enerd, wcycle);
1783 post_process_forces(cr, step, nrnb, wcycle,
1784 top, box, as_rvec_array(x.unpaddedArrayRef().data()), forceOut.f, &forceOut.forceWithVirial,
1785 vir_force, mdatoms, graph, fr, vsite,
1789 if (flags & GMX_FORCE_ENERGY)
1791 /* Sum the potential energy terms from group contributions */
1792 sum_epot(&(enerd->grpp), enerd->term);
1794 if (!EI_TPI(inputrec->eI))
1796 checkPotentialEnergyValidity(step, *enerd, *inputrec);
1802 void do_force(FILE *fplog,
1803 const t_commrec *cr,
1804 const gmx_multisim_t *ms,
1805 const t_inputrec *inputrec,
1807 gmx_enfrot *enforcedRotation,
1810 gmx_wallcycle_t wcycle,
1811 gmx_localtop_t *top,
1812 const gmx_groups_t *groups,
1814 gmx::ArrayRefWithPadding<gmx::RVec> x, //NOLINT(performance-unnecessary-value-param)
1816 gmx::ArrayRefWithPadding<gmx::RVec> force, //NOLINT(performance-unnecessary-value-param)
1818 const t_mdatoms *mdatoms,
1819 gmx_enerdata_t *enerd,
1821 gmx::ArrayRef<real> lambda,
1824 gmx::PpForceWorkload *ppForceWorkload,
1825 const gmx_vsite_t *vsite,
1830 const DDBalanceRegionHandler &ddBalanceRegionHandler)
1832 /* modify force flag if not doing nonbonded */
1833 if (!fr->bNonbonded)
1835 flags &= ~GMX_FORCE_NONBONDED;
1838 switch (inputrec->cutoff_scheme)
1841 do_force_cutsVERLET(fplog, cr, ms, inputrec,
1842 awh, enforcedRotation, step, nrnb, wcycle,
1849 lambda.data(), graph,
1856 ddBalanceRegionHandler);
1859 do_force_cutsGROUP(fplog, cr, ms, inputrec,
1860 awh, enforcedRotation, step, nrnb, wcycle,
1867 lambda.data(), graph,
1871 ddBalanceRegionHandler);
1874 gmx_incons("Invalid cut-off scheme passed!");
1877 /* In case we don't have constraints and are using GPUs, the next balancing
1878 * region starts here.
1879 * Some "special" work at the end of do_force_cuts?, such as vsite spread,
1880 * virial calculation and COM pulling, is not thus not included in
1881 * the balance timing, which is ok as most tasks do communication.
1883 ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);
1887 void do_constrain_first(FILE *fplog, gmx::Constraints *constr,
1888 const t_inputrec *ir, const t_mdatoms *md,
1891 int i, m, start, end;
1893 real dt = ir->delta_t;
1897 /* We need to allocate one element extra, since we might use
1898 * (unaligned) 4-wide SIMD loads to access rvec entries.
1900 snew(savex, state->natoms + 1);
1907 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
1908 start, md->homenr, end);
1910 /* Do a first constrain to reset particles... */
1911 step = ir->init_step;
1914 char buf[STEPSTRSIZE];
1915 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
1916 gmx_step_str(step, buf));
1920 /* constrain the current position */
1921 constr->apply(TRUE, FALSE,
1923 state->x.rvec_array(), state->x.rvec_array(), nullptr,
1925 state->lambda[efptBONDED], &dvdl_dum,
1926 nullptr, nullptr, gmx::ConstraintVariable::Positions);
1929 /* constrain the inital velocity, and save it */
1930 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1931 constr->apply(TRUE, FALSE,
1933 state->x.rvec_array(), state->v.rvec_array(), state->v.rvec_array(),
1935 state->lambda[efptBONDED], &dvdl_dum,
1936 nullptr, nullptr, gmx::ConstraintVariable::Velocities);
1938 /* constrain the inital velocities at t-dt/2 */
1939 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
1941 auto x = makeArrayRef(state->x).subArray(start, end);
1942 auto v = makeArrayRef(state->v).subArray(start, end);
1943 for (i = start; (i < end); i++)
1945 for (m = 0; (m < DIM); m++)
1947 /* Reverse the velocity */
1949 /* Store the position at t-dt in buf */
1950 savex[i][m] = x[i][m] + dt*v[i][m];
1953 /* Shake the positions at t=-dt with the positions at t=0
1954 * as reference coordinates.
1958 char buf[STEPSTRSIZE];
1959 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
1960 gmx_step_str(step, buf));
1963 constr->apply(TRUE, FALSE,
1965 state->x.rvec_array(), savex, nullptr,
1967 state->lambda[efptBONDED], &dvdl_dum,
1968 state->v.rvec_array(), nullptr, gmx::ConstraintVariable::Positions);
1970 for (i = start; i < end; i++)
1972 for (m = 0; m < DIM; m++)
1974 /* Re-reverse the velocities */
1982 void put_atoms_in_box_omp(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x)
1985 nth = gmx_omp_nthreads_get(emntDefault);
1987 #pragma omp parallel for num_threads(nth) schedule(static)
1988 for (t = 0; t < nth; t++)
1992 size_t natoms = x.size();
1993 size_t offset = (natoms*t )/nth;
1994 size_t len = (natoms*(t + 1))/nth - offset;
1995 put_atoms_in_box(ePBC, box, x.subArray(offset, len));
1997 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2001 void initialize_lambdas(FILE *fplog,
2002 const t_inputrec &ir,
2005 gmx::ArrayRef<real> lambda,
2008 /* TODO: Clean up initialization of fep_state and lambda in
2009 t_state. This function works, but could probably use a logic
2010 rewrite to keep all the different types of efep straight. */
2012 if ((ir.efep == efepNO) && (!ir.bSimTemp))
2017 const t_lambda *fep = ir.fepvals;
2020 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2021 if checkpoint is set -- a kludge is in for now
2025 for (int i = 0; i < efptNR; i++)
2028 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2029 if (fep->init_lambda >= 0) /* if it's -1, it was never initialized */
2031 thisLambda = fep->init_lambda;
2035 thisLambda = fep->all_lambda[i][fep->init_fep_state];
2039 lambda[i] = thisLambda;
2041 if (lam0 != nullptr)
2043 lam0[i] = thisLambda;
2048 /* need to rescale control temperatures to match current state */
2049 for (int i = 0; i < ir.opts.ngtc; i++)
2051 if (ir.opts.ref_t[i] > 0)
2053 ir.opts.ref_t[i] = ir.simtempvals->temperatures[fep->init_fep_state];
2058 /* Send to the log the information on the current lambdas */
2059 if (fplog != nullptr)
2061 fprintf(fplog, "Initial vector of lambda components:[ ");
2062 for (const auto &l : lambda)
2064 fprintf(fplog, "%10.4f ", l);
2066 fprintf(fplog, "]\n");