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) 2011-2019,2020,2021, 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.
39 * \brief Implements the integrator for normal molecular dynamics simulations
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
55 #include "gromacs/applied_forces/awh/awh.h"
56 #include "gromacs/applied_forces/awh/read_params.h"
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/collect.h"
59 #include "gromacs/domdec/dlbtiming.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_network.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localtopologychecker.h"
65 #include "gromacs/domdec/mdsetup.h"
66 #include "gromacs/domdec/partition.h"
67 #include "gromacs/essentialdynamics/edsam.h"
68 #include "gromacs/ewald/pme_load_balancing.h"
69 #include "gromacs/ewald/pme_pp.h"
70 #include "gromacs/fileio/trxio.h"
71 #include "gromacs/gmxlib/network.h"
72 #include "gromacs/gmxlib/nrnb.h"
73 #include "gromacs/gpu_utils/device_stream_manager.h"
74 #include "gromacs/gpu_utils/gpu_utils.h"
75 #include "gromacs/imd/imd.h"
76 #include "gromacs/listed_forces/listed_forces.h"
77 #include "gromacs/math/functions.h"
78 #include "gromacs/math/invertmatrix.h"
79 #include "gromacs/math/vec.h"
80 #include "gromacs/math/vectypes.h"
81 #include "gromacs/mdlib/checkpointhandler.h"
82 #include "gromacs/mdlib/compute_io.h"
83 #include "gromacs/mdlib/constr.h"
84 #include "gromacs/mdlib/coupling.h"
85 #include "gromacs/mdlib/ebin.h"
86 #include "gromacs/mdlib/enerdata_utils.h"
87 #include "gromacs/mdlib/energyoutput.h"
88 #include "gromacs/mdlib/expanded.h"
89 #include "gromacs/mdlib/force.h"
90 #include "gromacs/mdlib/force_flags.h"
91 #include "gromacs/mdlib/forcerec.h"
92 #include "gromacs/mdlib/freeenergyparameters.h"
93 #include "gromacs/mdlib/md_support.h"
94 #include "gromacs/mdlib/mdatoms.h"
95 #include "gromacs/mdlib/mdoutf.h"
96 #include "gromacs/mdlib/membed.h"
97 #include "gromacs/mdlib/resethandler.h"
98 #include "gromacs/mdlib/sighandler.h"
99 #include "gromacs/mdlib/simulationsignal.h"
100 #include "gromacs/mdlib/stat.h"
101 #include "gromacs/mdlib/stophandler.h"
102 #include "gromacs/mdlib/tgroup.h"
103 #include "gromacs/mdlib/trajectory_writing.h"
104 #include "gromacs/mdlib/update.h"
105 #include "gromacs/mdlib/update_constrain_gpu.h"
106 #include "gromacs/mdlib/update_vv.h"
107 #include "gromacs/mdlib/vcm.h"
108 #include "gromacs/mdlib/vsite.h"
109 #include "gromacs/mdrunutility/freeenergy.h"
110 #include "gromacs/mdrunutility/handlerestart.h"
111 #include "gromacs/mdrunutility/multisim.h"
112 #include "gromacs/mdrunutility/printtime.h"
113 #include "gromacs/mdtypes/awh_history.h"
114 #include "gromacs/mdtypes/awh_params.h"
115 #include "gromacs/mdtypes/commrec.h"
116 #include "gromacs/mdtypes/df_history.h"
117 #include "gromacs/mdtypes/energyhistory.h"
118 #include "gromacs/mdtypes/fcdata.h"
119 #include "gromacs/mdtypes/forcebuffers.h"
120 #include "gromacs/mdtypes/forcerec.h"
121 #include "gromacs/mdtypes/group.h"
122 #include "gromacs/mdtypes/inputrec.h"
123 #include "gromacs/mdtypes/interaction_const.h"
124 #include "gromacs/mdtypes/md_enums.h"
125 #include "gromacs/mdtypes/mdatom.h"
126 #include "gromacs/mdtypes/mdrunoptions.h"
127 #include "gromacs/mdtypes/multipletimestepping.h"
128 #include "gromacs/mdtypes/observableshistory.h"
129 #include "gromacs/mdtypes/pullhistory.h"
130 #include "gromacs/mdtypes/simulation_workload.h"
131 #include "gromacs/mdtypes/state.h"
132 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
133 #include "gromacs/modularsimulator/energydata.h"
134 #include "gromacs/nbnxm/gpu_data_mgmt.h"
135 #include "gromacs/nbnxm/nbnxm.h"
136 #include "gromacs/pbcutil/pbc.h"
137 #include "gromacs/pulling/output.h"
138 #include "gromacs/pulling/pull.h"
139 #include "gromacs/swap/swapcoords.h"
140 #include "gromacs/timing/wallcycle.h"
141 #include "gromacs/timing/walltime_accounting.h"
142 #include "gromacs/topology/atoms.h"
143 #include "gromacs/topology/idef.h"
144 #include "gromacs/topology/mtop_util.h"
145 #include "gromacs/topology/topology.h"
146 #include "gromacs/trajectory/trajectoryframe.h"
147 #include "gromacs/utility/basedefinitions.h"
148 #include "gromacs/utility/cstringutil.h"
149 #include "gromacs/utility/fatalerror.h"
150 #include "gromacs/utility/logger.h"
151 #include "gromacs/utility/real.h"
152 #include "gromacs/utility/smalloc.h"
154 #include "legacysimulator.h"
155 #include "replicaexchange.h"
158 using gmx::SimulationSignaller;
160 void gmx::LegacySimulator::do_md()
162 // TODO Historically, the EM and MD "integrators" used different
163 // names for the t_inputrec *parameter, but these must have the
164 // same name, now that it's a member of a struct. We use this ir
165 // alias to avoid a large ripple of nearly useless changes.
166 // t_inputrec is being replaced by IMdpOptionsProvider, so this
167 // will go away eventually.
168 const t_inputrec* ir = inputrec;
170 int64_t step, step_rel;
171 double t, t0 = ir->init_t;
172 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
173 gmx_bool bNS = FALSE, bNStList, bStopCM, bFirstStep, bInitStep, bLastStep = FALSE;
174 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
175 gmx_bool do_ene, do_log, do_verbose;
176 gmx_bool bMasterState;
177 unsigned int force_flags;
178 tensor force_vir = { { 0 } }, shake_vir = { { 0 } }, total_vir = { { 0 } }, pres = { { 0 } };
181 matrix pressureCouplingMu, M;
182 gmx_repl_ex_t repl_ex = nullptr;
183 gmx_global_stat_t gstat;
184 gmx_shellfc_t* shellfc;
185 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
188 std::vector<RVec> cbuf;
193 real saved_conserved_quantity = 0;
196 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
198 /* PME load balancing data for GPU kernels */
199 gmx_bool bPMETune = FALSE;
200 gmx_bool bPMETunePrinting = FALSE;
202 bool bInteractiveMDstep = false;
204 SimulationSignals signals;
205 // Most global communnication stages don't propagate mdrun
206 // signals, and will use this object to achieve that.
207 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
209 if (!mdrunOptions.writeConfout)
211 // This is on by default, and the main known use case for
212 // turning it off is for convenience in benchmarking, which is
213 // something that should not show up in the general user
218 "The -noconfout functionality is deprecated, and may be removed in a "
222 /* md-vv uses averaged full step velocities for T-control
223 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
224 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
225 bTrotter = (EI_VV(ir->eI)
226 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
228 const bool bRerunMD = false;
230 int nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
231 bGStatEveryStep = (nstglobalcomm == 1);
233 const SimulationGroups* groups = &top_global.groups;
235 std::unique_ptr<EssentialDynamics> ed = nullptr;
236 if (opt2bSet("-ei", nfile, fnm))
238 /* Initialize essential dynamics sampling */
239 ed = init_edsam(mdlog,
240 opt2fn_null("-ei", nfile, fnm),
241 opt2fn("-eo", nfile, fnm),
251 else if (observablesHistory->edsamHistory)
254 "The checkpoint is from a run with essential dynamics sampling, "
255 "but the current run did not specify the -ei option. "
256 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
259 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
260 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
261 initialize_lambdas(fplog,
265 ir->simtempvals->temperatures,
266 gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
270 Update upd(*ir, deform);
271 bool doSimulatedAnnealing = false;
273 // TODO: Avoid changing inputrec (#3854)
274 // Simulated annealing updates the reference temperature.
275 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
276 doSimulatedAnnealing = initSimulatedAnnealing(nonConstInputrec, &upd);
278 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
280 const t_fcdata& fcdata = *fr->fcdata;
282 bool simulationsShareState = false;
283 int nstSignalComm = nstglobalcomm;
285 // TODO This implementation of ensemble orientation restraints is nasty because
286 // a user can't just do multi-sim with single-sim orientation restraints.
287 bool usingEnsembleRestraints = (fcdata.disres->nsystems > 1) || ((ms != nullptr) && fcdata.orires);
288 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim() && (ms != nullptr));
290 // Replica exchange, ensemble restraints and AWH need all
291 // simulations to remain synchronized, so they need
292 // checkpoints and stop conditions to act on the same step, so
293 // the propagation of such signals must take place between
294 // simulations, not just within simulations.
295 // TODO: Make algorithm initializers set these flags.
296 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
298 if (simulationsShareState)
300 // Inter-simulation signal communication does not need to happen
301 // often, so we use a minimum of 200 steps to reduce overhead.
302 const int c_minimumInterSimulationSignallingInterval = 200;
303 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
308 if (startingBehavior != StartingBehavior::RestartWithAppending)
310 pleaseCiteCouplingAlgorithms(fplog, *ir);
312 gmx_mdoutf* outf = init_mdoutf(fplog,
324 simulationsShareState,
326 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
330 mdoutf_get_fp_dhdl(outf),
333 simulationsShareState,
336 gstat = global_stat_init(ir);
338 const auto& simulationWork = runScheduleWork->simulationWork;
339 const bool useGpuForPme = simulationWork.useGpuPme;
340 const bool useGpuForNonbonded = simulationWork.useGpuNonbonded;
341 const bool useGpuForBufferOps = simulationWork.useGpuBufferOps;
342 const bool useGpuForUpdate = simulationWork.useGpuUpdate;
344 /* Check for polarizable models and flexible constraints */
345 shellfc = init_shell_flexcon(fplog,
347 constr ? constr->numFlexibleConstraints() : 0,
353 double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
354 if ((io > 2000) && MASTER(cr))
356 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
360 // Local state only becomes valid now.
361 std::unique_ptr<t_state> stateInstance;
364 gmx_localtop_t top(top_global.ffparams);
366 ForceBuffers f(simulationWork.useMts,
367 ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
368 ? PinningPolicy::PinnedIfSupported
369 : PinningPolicy::CannotBePinned);
370 const t_mdatoms* md = mdAtoms->mdatoms();
371 if (DOMAINDECOMP(cr))
373 stateInstance = std::make_unique<t_state>();
374 state = stateInstance.get();
375 dd_init_local_state(*cr->dd, state_global, state);
377 /* Distribute the charge groups over the nodes from the master node */
378 dd_partition_system(fplog,
399 upd.updateAfterPartition(state->natoms,
400 md->cFREEZE ? gmx::arrayRefFromArray(md->cFREEZE, md->nr)
401 : gmx::ArrayRef<const unsigned short>(),
402 md->cTC ? gmx::arrayRefFromArray(md->cTC, md->nr)
403 : gmx::ArrayRef<const unsigned short>());
407 state_change_natoms(state_global, state_global->natoms);
408 /* Copy the pointer to the global state */
409 state = state_global;
411 /* Generate and initialize new topology */
412 mdAlgorithmsSetupAtomData(cr, *ir, top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
414 upd.updateAfterPartition(state->natoms,
415 md->cFREEZE ? gmx::arrayRefFromArray(md->cFREEZE, md->nr)
416 : gmx::ArrayRef<const unsigned short>(),
417 md->cTC ? gmx::arrayRefFromArray(md->cTC, md->nr)
418 : gmx::ArrayRef<const unsigned short>());
421 std::unique_ptr<UpdateConstrainGpu> integrator;
423 StatePropagatorDataGpu* stateGpu = fr->stateGpu;
425 // TODO: the assertions below should be handled by UpdateConstraintsBuilder.
428 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
429 || constr->numConstraintsTotal() == 0,
430 "Constraints in domain decomposition are only supported with update "
431 "groups if using GPU update.\n");
432 GMX_RELEASE_ASSERT(ir->eConstrAlg != ConstraintAlgorithm::Shake || constr == nullptr
433 || constr->numConstraintsTotal() == 0,
434 "SHAKE is not supported with GPU update.");
435 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
436 "Either PME or short-ranged non-bonded interaction tasks must run on "
437 "the GPU to use GPU update.\n");
438 GMX_RELEASE_ASSERT(ir->eI == IntegrationAlgorithm::MD,
439 "Only the md integrator is supported with the GPU update.\n");
441 ir->etc != TemperatureCoupling::NoseHoover,
442 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
444 ir->epc == PressureCoupling::No || ir->epc == PressureCoupling::ParrinelloRahman
445 || ir->epc == PressureCoupling::Berendsen || ir->epc == PressureCoupling::CRescale,
446 "Only Parrinello-Rahman, Berendsen, and C-rescale pressure coupling are supported "
447 "with the GPU update.\n");
448 GMX_RELEASE_ASSERT(!md->haveVsites,
449 "Virtual sites are not supported with the GPU update.\n");
450 GMX_RELEASE_ASSERT(ed == nullptr,
451 "Essential dynamics is not supported with the GPU update.\n");
452 GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(*ir->pull),
453 "Constraints pulling is not supported with the GPU update.\n");
454 GMX_RELEASE_ASSERT(fcdata.orires == nullptr,
455 "Orientation restraints are not supported with the GPU update.\n");
457 ir->efep == FreeEnergyPerturbationType::No
458 || (!haveFepPerturbedMasses(top_global) && !havePerturbedConstraints(top_global)),
459 "Free energy perturbation of masses and constraints are not supported with the GPU "
462 if (constr != nullptr && constr->numConstraintsTotal() > 0)
466 .appendText("Updating coordinates and applying constraints on the GPU.");
470 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
472 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
473 "Device stream manager should be initialized in order to use GPU "
474 "update-constraints.");
476 fr->deviceStreamManager->streamIsValid(gmx::DeviceStreamType::UpdateAndConstraints),
477 "Update stream should be initialized in order to use GPU "
478 "update-constraints.");
479 integrator = std::make_unique<UpdateConstrainGpu>(
483 fr->deviceStreamManager->context(),
484 fr->deviceStreamManager->stream(gmx::DeviceStreamType::UpdateAndConstraints),
485 stateGpu->xUpdatedOnDevice(),
488 integrator->setPbc(PbcType::Xyz, state->box);
491 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
493 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
497 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
500 // NOTE: The global state is no longer used at this point.
501 // But state_global is still used as temporary storage space for writing
502 // the global state to file and potentially for replica exchange.
503 // (Global topology should persist.)
505 update_mdatoms(mdAtoms->mdatoms(), state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
509 /* Check nstexpanded here, because the grompp check was broken */
510 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
513 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
515 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
520 EnergyData::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
523 preparePrevStepPullCom(ir,
525 gmx::arrayRefFromArray(md->massT, md->nr),
529 startingBehavior != StartingBehavior::NewSimulation);
531 // TODO: Remove this by converting AWH into a ForceProvider
532 auto awh = prepareAwhModule(fplog,
537 startingBehavior != StartingBehavior::NewSimulation,
539 opt2fn("-awh", nfile, fnm),
542 if (useReplicaExchange && MASTER(cr))
544 repl_ex = init_replica_exchange(fplog, ms, top_global.natoms, ir, replExParams);
546 /* PME tuning is only supported in the Verlet scheme, with PME for
547 * Coulomb. It is not supported with only LJ PME. */
548 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
549 && ir->cutoff_scheme != CutoffScheme::Group);
551 pme_load_balancing_t* pme_loadbal = nullptr;
555 &pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata, fr->nbv->useGpu());
558 if (!ir->bContinuation)
560 if (state->flags & enumValueToBitMask(StateEntry::V))
562 auto v = makeArrayRef(state->v);
563 /* Set the velocities of vsites, shells and frozen atoms to zero */
564 for (i = 0; i < md->homenr; i++)
566 if (md->ptype[i] == ParticleType::Shell)
570 else if (md->cFREEZE)
572 for (m = 0; m < DIM; m++)
574 if (ir->opts.nFreeze[md->cFREEZE[i]][m])
585 /* Constrain the initial coordinates and velocities */
586 do_constrain_first(fplog,
591 state->x.arrayRefWithPadding(),
592 state->v.arrayRefWithPadding(),
594 state->lambda[FreeEnergyPerturbationCouplingType::Bonded]);
598 const int nstfep = computeFepPeriod(*ir, replExParams);
600 /* Be REALLY careful about what flags you set here. You CANNOT assume
601 * this is the first step, since we might be restarting from a checkpoint,
602 * and in that case we should not do any modifications to the state.
604 bStopCM = (ir->comm_mode != ComRemovalAlgorithm::No && !ir->bContinuation);
606 // When restarting from a checkpoint, it can be appropriate to
607 // initialize ekind from quantities in the checkpoint. Otherwise,
608 // compute_globals must initialize ekind before the simulation
609 // starts/restarts. However, only the master rank knows what was
610 // found in the checkpoint file, so we have to communicate in
611 // order to coordinate the restart.
613 // TODO Consider removing this communication if/when checkpoint
614 // reading directly follows .tpr reading, because all ranks can
615 // agree on hasReadEkinState at that time.
616 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
619 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr->mpi_comm_mygroup);
621 if (hasReadEkinState)
623 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
626 unsigned int cglo_flags =
627 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
628 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
630 bSumEkinhOld = FALSE;
632 t_vcm vcm(top_global.groups, *ir);
633 reportComRemovalInfo(fplog, vcm);
635 /* To minimize communication, compute_globals computes the COM velocity
636 * and the kinetic energy for the velocities without COM motion removed.
637 * Thus to get the kinetic energy without the COM contribution, we need
638 * to call compute_globals twice.
640 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
642 unsigned int cglo_flags_iteration = cglo_flags;
643 if (bStopCM && cgloIteration == 0)
645 cglo_flags_iteration |= CGLO_STOPCM;
646 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
648 if (DOMAINDECOMP(cr) && dd_localTopologyChecker(*cr->dd).shouldCheckNumberOfBondedInteractions()
649 && cgloIteration == 0)
651 cglo_flags_iteration |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
653 compute_globals(gstat,
658 makeConstArrayRef(state->x),
659 makeConstArrayRef(state->v),
670 gmx::ArrayRef<real>{},
674 cglo_flags_iteration);
675 if (cglo_flags_iteration & CGLO_STOPCM)
677 /* At initialization, do not pass x with acceleration-correction mode
678 * to avoid (incorrect) correction of the initial coordinates.
680 auto x = (vcm.mode == ComRemovalAlgorithm::LinearAccelerationCorrection)
682 : makeArrayRef(state->x);
683 process_and_stopcm_grp(fplog, &vcm, *md, x, makeArrayRef(state->v));
684 inc_nrnb(nrnb, eNR_STOPCM, md->homenr);
687 if (DOMAINDECOMP(cr))
689 dd_localTopologyChecker(cr->dd)->checkNumberOfBondedInteractions(
690 &top, makeConstArrayRef(state->x), state->box);
692 if (ir->eI == IntegrationAlgorithm::VVAK)
694 /* a second call to get the half step temperature initialized as well */
695 /* we do the same call as above, but turn the pressure off -- internally to
696 compute_globals, this is recognized as a velocity verlet half-step
697 kinetic energy calculation. This minimized excess variables, but
698 perhaps loses some logic?*/
700 compute_globals(gstat,
705 makeConstArrayRef(state->x),
706 makeConstArrayRef(state->v),
717 gmx::ArrayRef<real>{},
721 cglo_flags & ~CGLO_PRESSURE);
724 /* Calculate the initial half step temperature, and save the ekinh_old */
725 if (startingBehavior == StartingBehavior::NewSimulation)
727 for (i = 0; (i < ir->opts.ngtc); i++)
729 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
733 /* need to make an initiation call to get the Trotter variables set, as well as other constants
734 for non-trotter temperature control */
735 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
739 if (!ir->bContinuation)
741 if (constr && ir->eConstrAlg == ConstraintAlgorithm::Lincs)
744 "RMS relative constraint deviation after constraining: %.2e\n",
747 if (EI_STATE_VELOCITY(ir->eI))
749 real temp = enerd->term[F_TEMP];
750 if (ir->eI != IntegrationAlgorithm::VV)
752 /* Result of Ekin averaged over velocities of -half
753 * and +half step, while we only have -half step here.
757 fprintf(fplog, "Initial temperature: %g K\n", temp);
762 fprintf(stderr, "starting mdrun '%s'\n", *(top_global.name));
765 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
769 sprintf(tbuf, "%s", "infinite");
771 if (ir->init_step > 0)
774 "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
775 gmx_step_str(ir->init_step + ir->nsteps, sbuf),
777 gmx_step_str(ir->init_step, sbuf2),
778 ir->init_step * ir->delta_t);
782 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
784 fprintf(fplog, "\n");
787 walltime_accounting_start_time(walltime_accounting);
788 wallcycle_start(wcycle, WallCycleCounter::Run);
789 print_start(fplog, cr, walltime_accounting, "mdrun");
791 /***********************************************************
795 ************************************************************/
798 /* Skip the first Nose-Hoover integration when we get the state from tpx */
799 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
800 bSumEkinhOld = FALSE;
802 bNeedRepartition = FALSE;
804 step = ir->init_step;
807 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
808 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
809 simulationsShareState,
812 mdrunOptions.reproducible,
814 mdrunOptions.maximumHoursToRun,
819 walltime_accounting);
821 auto checkpointHandler = std::make_unique<CheckpointHandler>(
822 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]),
823 simulationsShareState,
826 mdrunOptions.writeConfout,
827 mdrunOptions.checkpointOptions.period);
829 const bool resetCountersIsLocal = true;
830 auto resetHandler = std::make_unique<ResetHandler>(
831 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
832 !resetCountersIsLocal,
835 mdrunOptions.timingOptions.resetHalfway,
836 mdrunOptions.maximumHoursToRun,
839 walltime_accounting);
841 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
843 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
845 logInitialMultisimStatus(ms, cr, mdlog, simulationsShareState, ir->nsteps, ir->init_step);
848 /* and stop now if we should */
849 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
853 /* Determine if this is a neighbor search step */
854 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
856 if (bPMETune && bNStList)
858 // This has to be here because PME load balancing is called so early.
859 // TODO: Move to after all booleans are defined.
860 if (useGpuForUpdate && !bFirstStep)
862 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
863 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
865 /* PME grid + cut-off optimization with GPUs or PME nodes */
866 pme_loadbal_do(pme_loadbal,
868 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
879 simulationWork.useGpuPmePpCommunication);
882 wallcycle_start(wcycle, WallCycleCounter::Step);
884 bLastStep = (step_rel == ir->nsteps);
885 t = t0 + step * ir->delta_t;
887 // TODO Refactor this, so that nstfep does not need a default value of zero
888 if (ir->efep != FreeEnergyPerturbationType::No || ir->bSimTemp)
890 /* find and set the current lambdas */
891 state->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
893 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
894 bDoFEP = ((ir->efep != FreeEnergyPerturbationType::No) && do_per_step(step, nstfep));
895 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
899 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
900 && do_per_step(step, replExParams.exchangeInterval));
902 if (doSimulatedAnnealing)
904 // TODO: Avoid changing inputrec (#3854)
905 // Simulated annealing updates the reference temperature.
906 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
907 update_annealing_target_temp(nonConstInputrec, t, &upd);
910 /* Stop Center of Mass motion */
911 bStopCM = (ir->comm_mode != ComRemovalAlgorithm::No && do_per_step(step, ir->nstcomm));
913 /* Determine whether or not to do Neighbour Searching */
914 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
916 /* Note that the stopHandler will cause termination at nstglobalcomm
917 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
918 * nstpcouple steps, we have computed the half-step kinetic energy
919 * of the previous step and can always output energies at the last step.
921 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
923 /* do_log triggers energy and virial calculation. Because this leads
924 * to different code paths, forces can be different. Thus for exact
925 * continuation we should avoid extra log output.
926 * Note that the || bLastStep can result in non-exact continuation
927 * beyond the last step. But we don't consider that to be an issue.
929 do_log = (do_per_step(step, ir->nstlog)
930 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
931 do_verbose = mdrunOptions.verbose
932 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
934 if (useGpuForUpdate && !bFirstStep && bNS)
936 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
937 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
938 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
939 // Copy coordinate from the GPU when needed at the search step.
940 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
941 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
942 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
943 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
946 // We only need to calculate virtual velocities if we are writing them in the current step
947 const bool needVirtualVelocitiesThisStep =
949 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep());
951 if (vsite != nullptr)
953 // Virtual sites need to be updated before domain decomposition and forces are calculated
954 wallcycle_start(wcycle, WallCycleCounter::VsiteConstr);
955 // md-vv calculates virtual velocities once it has full-step real velocities
956 vsite->construct(state->x,
959 (!EI_VV(inputrec->eI) && needVirtualVelocitiesThisStep)
960 ? VSiteOperation::PositionsAndVelocities
961 : VSiteOperation::Positions);
962 wallcycle_stop(wcycle, WallCycleCounter::VsiteConstr);
965 if (bNS && !(bFirstStep && ir->bContinuation))
967 bMasterState = FALSE;
968 /* Correct the new box if it is too skewed */
969 if (inputrecDynamicBox(ir))
971 if (correct_box(fplog, step, state->box))
974 // If update is offloaded, it should be informed about the box size change
977 integrator->setPbc(PbcType::Xyz, state->box);
981 if (DOMAINDECOMP(cr) && bMasterState)
983 dd_collect_state(cr->dd, state, state_global);
986 if (DOMAINDECOMP(cr))
988 /* Repartition the domain decomposition */
989 dd_partition_system(fplog,
1009 do_verbose && !bPMETunePrinting);
1010 upd.updateAfterPartition(state->natoms,
1011 md->cFREEZE ? gmx::arrayRefFromArray(md->cFREEZE, md->nr)
1012 : gmx::ArrayRef<const unsigned short>(),
1013 md->cTC ? gmx::arrayRefFromArray(md->cTC, md->nr)
1014 : gmx::ArrayRef<const unsigned short>());
1018 // Allocate or re-size GPU halo exchange object, if necessary
1019 if (bNS && havePPDomainDecomposition(cr) && simulationWork.useGpuHaloExchange)
1021 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
1022 "GPU device manager has to be initialized to use GPU "
1023 "version of halo exchange.");
1024 constructGpuHaloExchange(mdlog, *cr, *fr->deviceStreamManager, wcycle);
1027 if (MASTER(cr) && do_log)
1029 gmx::EnergyOutput::printHeader(
1030 fplog, step, t); /* can we improve the information printed here? */
1033 if (ir->efep != FreeEnergyPerturbationType::No)
1035 update_mdatoms(mdAtoms->mdatoms(), state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
1040 /* We need the kinetic energy at minus the half step for determining
1041 * the full step kinetic energy and possibly for T-coupling.*/
1042 /* This may not be quite working correctly yet . . . . */
1043 int cglo_flags = CGLO_GSTAT | CGLO_TEMPERATURE;
1044 if (DOMAINDECOMP(cr) && dd_localTopologyChecker(*cr->dd).shouldCheckNumberOfBondedInteractions())
1046 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
1048 compute_globals(gstat,
1053 makeConstArrayRef(state->x),
1054 makeConstArrayRef(state->v),
1065 gmx::ArrayRef<real>{},
1070 if (DOMAINDECOMP(cr))
1072 dd_localTopologyChecker(cr->dd)->checkNumberOfBondedInteractions(
1073 &top, makeConstArrayRef(state->x), state->box);
1076 clear_mat(force_vir);
1078 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
1080 /* Determine the energy and pressure:
1081 * at nstcalcenergy steps and at energy output steps (set below).
1083 if (EI_VV(ir->eI) && (!bInitStep))
1085 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1086 bCalcVir = bCalcEnerStep
1087 || (ir->epc != PressureCoupling::No
1088 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
1092 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1093 bCalcVir = bCalcEnerStep
1094 || (ir->epc != PressureCoupling::No && do_per_step(step, ir->nstpcouple));
1096 bCalcEner = bCalcEnerStep;
1098 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1100 if (do_ene || do_log || bDoReplEx)
1106 /* Do we need global communication ? */
1107 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
1108 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
1110 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
1111 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
1112 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
1113 if (simulationWork.useMts && !do_per_step(step, ir->nstfout))
1115 // TODO: merge this with stepWork.useOnlyMtsCombinedForceBuffer
1116 force_flags |= GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE;
1121 /* Now is the time to relax the shells */
1122 relax_shell_flexcon(fplog,
1125 mdrunOptions.verbose,
1137 state->x.arrayRefWithPadding(),
1138 state->v.arrayRefWithPadding(),
1153 ddBalanceRegionHandler);
1157 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
1158 is updated (or the AWH update will be performed twice for one step when continuing).
1159 It would be best to call this update function from do_md_trajectory_writing but that
1160 would occur after do_force. One would have to divide the update_awh function into one
1161 function applying the AWH force and one doing the AWH bias update. The update AWH
1162 bias function could then be called after do_md_trajectory_writing (then containing
1163 update_awh_history). The checkpointing will in the future probably moved to the start
1164 of the md loop which will rid of this issue. */
1165 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
1167 awh->updateHistory(state_global->awhHistory.get());
1170 /* The coordinates (x) are shifted (to get whole molecules)
1172 * This is parallellized as well, and does communication too.
1173 * Check comments in sim_util.c
1188 state->x.arrayRefWithPadding(),
1200 ed ? ed->getLegacyED() : nullptr,
1201 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1202 ddBalanceRegionHandler);
1205 // VV integrators do not need the following velocity half step
1206 // if it is the first step after starting from a checkpoint.
1207 // That is, the half step is needed on all other steps, and
1208 // also the first step when starting from a .tpr file.
1211 integrateVVFirstStep(step,
1243 &saved_conserved_quantity,
1252 if (vsite != nullptr && needVirtualVelocitiesThisStep)
1254 // Positions were calculated earlier
1255 wallcycle_start(wcycle, WallCycleCounter::VsiteConstr);
1256 vsite->construct(state->x, state->v, state->box, VSiteOperation::Velocities);
1257 wallcycle_stop(wcycle, WallCycleCounter::VsiteConstr);
1261 /* ######## END FIRST UPDATE STEP ############## */
1262 /* ######## If doing VV, we now have v(dt) ###### */
1265 /* perform extended ensemble sampling in lambda - we don't
1266 actually move to the new state before outputting
1267 statistics, but if performing simulated tempering, we
1268 do update the velocities and the tau_t. */
1269 // TODO: Avoid changing inputrec (#3854)
1270 // Simulated tempering updates the reference temperature.
1271 // Expanded ensemble without simulated tempering does not change the inputrec.
1272 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
1273 lamnew = ExpandedEnsembleDynamics(fplog,
1281 state->v.rvec_array(),
1283 md->cTC ? gmx::arrayRefFromArray(md->cTC, md->nr)
1284 : gmx::ArrayRef<const unsigned short>());
1285 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1288 copy_df_history(state_global->dfhist, state->dfhist);
1292 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1293 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1294 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1295 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1296 || checkpointHandler->isCheckpointingStep()))
1298 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1299 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1301 // Copy velocities if needed for the output/checkpointing.
1302 // NOTE: Copy on the search steps is done at the beginning of the step.
1303 if (useGpuForUpdate && !bNS
1304 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1306 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1307 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1309 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1310 // and update is offloaded hence forces are kept on the GPU for update and have not been
1311 // already transferred in do_force().
1312 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1313 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1314 // prior to GPU update.
1315 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1316 // copy call in do_force(...).
1317 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1318 // on host after the D2H copy in do_force(...).
1319 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1320 && do_per_step(step, ir->nstfout))
1322 stateGpu->copyForcesFromGpu(f.view().force(), AtomLocality::Local);
1323 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1325 /* Now we have the energies and forces corresponding to the
1326 * coordinates at time t. We must output all of this before
1329 do_md_trajectory_writing(fplog,
1346 checkpointHandler->isCheckpointingStep(),
1349 mdrunOptions.writeConfout,
1351 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1352 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x, t);
1354 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1355 if (startingBehavior != StartingBehavior::NewSimulation && bFirstStep
1356 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1358 copy_mat(state->svir_prev, shake_vir);
1359 copy_mat(state->fvir_prev, force_vir);
1362 stopHandler->setSignal();
1363 resetHandler->setSignal(walltime_accounting);
1365 if (bGStat || !PAR(cr))
1367 /* In parallel we only have to check for checkpointing in steps
1368 * where we do global communication,
1369 * otherwise the other nodes don't know.
1371 checkpointHandler->setSignal(walltime_accounting);
1374 /* ######### START SECOND UPDATE STEP ################# */
1376 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1377 controlled in preprocessing */
1379 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1381 gmx_bool bIfRandomize;
1382 bIfRandomize = update_randomize_velocities(ir,
1386 md->cTC ? gmx::arrayRefFromArray(md->cTC, md->nr)
1387 : gmx::ArrayRef<const unsigned short>(),
1388 gmx::arrayRefFromArray(md->invmass, md->nr),
1392 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1393 if (constr && bIfRandomize)
1395 constrain_velocities(constr, do_log, do_ene, step, state, nullptr, false, nullptr);
1398 /* Box is changed in update() when we do pressure coupling,
1399 * but we should still use the old box for energy corrections and when
1400 * writing it to the energy file, so it matches the trajectory files for
1401 * the same timestep above. Make a copy in a separate array.
1403 copy_mat(state->box, lastbox);
1407 if (!useGpuForUpdate)
1409 wallcycle_start(wcycle, WallCycleCounter::Update);
1411 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1421 md->cTC ? gmx::arrayRefFromArray(md->cTC, md->nr)
1422 : gmx::ArrayRef<const unsigned short>(),
1423 gmx::arrayRefFromArray(md->invmass, md->nr),
1426 TrotterSequence::Three);
1427 /* We can only do Berendsen coupling after we have summed
1428 * the kinetic energy or virial. Since the happens
1429 * in global_state after update, we should only do it at
1430 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1435 update_tcouple(step,
1441 md->cTC ? gmx::arrayRefFromArray(md->cTC, md->nr)
1442 : gmx::ArrayRef<const unsigned short>());
1443 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1446 /* With leap-frog type integrators we compute the kinetic energy
1447 * at a whole time step as the average of the half-time step kinetic
1448 * energies of two subsequent steps. Therefore we need to compute the
1449 * half step kinetic energy also if we need energies at the next step.
1451 const bool needHalfStepKineticEnergy =
1452 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1454 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1455 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1456 const bool doParrinelloRahman = (ir->epc == PressureCoupling::ParrinelloRahman
1457 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1461 GMX_ASSERT(!useGpuForUpdate, "GPU update is not supported with VVAK integrator.");
1463 integrateVVSecondStep(step,
1499 if (useGpuForUpdate)
1501 if (bNS && (bFirstStep || DOMAINDECOMP(cr)))
1503 integrator->set(stateGpu->getCoordinates(),
1504 stateGpu->getVelocities(),
1505 stateGpu->getForces(),
1509 // Copy data to the GPU after buffers might have being reinitialized
1510 /* The velocity copy is redundant if we had Center-of-Mass motion removed on
1511 * the previous step. We don't check that now. */
1512 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1513 // coordinates have been copied already if PME or buffer ops has not needed it this step.
1514 const bool useGpuPmeOnThisRank = runScheduleWork->simulationWork.useGpuPme
1515 && thisRankHasDuty(cr, DUTY_PME)
1516 && runScheduleWork->stepWork.computeSlowForces;
1517 if (!useGpuPmeOnThisRank && !runScheduleWork->stepWork.useGpuXBufferOps)
1519 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1523 if (simulationWork.useGpuPme && !runScheduleWork->simulationWork.useGpuPmePpCommunication
1524 && !thisRankHasDuty(cr, DUTY_PME))
1526 // The PME forces were recieved to the host, so have to be copied
1527 stateGpu->copyForcesToGpu(f.view().force(), AtomLocality::All);
1529 else if (!runScheduleWork->stepWork.useGpuFBufferOps)
1531 // The buffer ops were not offloaded this step, so the forces are on the
1532 // host and have to be copied
1533 stateGpu->copyForcesToGpu(f.view().force(), AtomLocality::Local);
1536 const bool doTemperatureScaling =
1537 (ir->etc != TemperatureCoupling::No
1538 && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1540 // This applies Leap-Frog, LINCS and SETTLE in succession
1541 integrator->integrate(
1542 stateGpu->getForcesReadyOnDeviceEvent(
1543 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1548 doTemperatureScaling,
1551 ir->nstpcouple * ir->delta_t,
1554 // Copy velocities D2H after update if:
1555 // - Globals are computed this step (includes the energy output steps).
1556 // - Temperature is needed for the next step.
1557 if (bGStat || needHalfStepKineticEnergy)
1559 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1560 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1565 /* With multiple time stepping we need to do an additional normal
1566 * update step to obtain the virial, as the actual MTS integration
1567 * using an acceleration where the slow forces are multiplied by mtsFactor.
1568 * Using that acceleration would result in a virial with the slow
1569 * force contribution would be a factor mtsFactor too large.
1571 if (simulationWork.useMts && bCalcVir && constr != nullptr)
1573 upd.update_for_constraint_virial(*ir,
1575 md->havePartiallyFrozenAtoms,
1576 gmx::arrayRefFromArray(md->invmass, md->nr),
1577 gmx::arrayRefFromArray(md->invMassPerDim, md->nr),
1579 f.view().forceWithPadding(),
1582 constrain_coordinates(constr,
1587 upd.xp()->arrayRefWithPadding(),
1593 ArrayRefWithPadding<const RVec> forceCombined =
1594 (simulationWork.useMts && step % ir->mtsLevels[1].stepFactor == 0)
1595 ? f.view().forceMtsCombinedWithPadding()
1596 : f.view().forceWithPadding();
1597 upd.update_coords(*ir,
1600 md->havePartiallyFrozenAtoms,
1601 gmx::arrayRefFromArray(md->ptype, md->nr),
1602 gmx::arrayRefFromArray(md->invmass, md->nr),
1603 gmx::arrayRefFromArray(md->invMassPerDim, md->nr),
1613 wallcycle_stop(wcycle, WallCycleCounter::Update);
1615 constrain_coordinates(constr,
1620 upd.xp()->arrayRefWithPadding(),
1622 bCalcVir && !simulationWork.useMts,
1625 upd.update_sd_second_half(*ir,
1629 gmx::arrayRefFromArray(md->ptype, md->nr),
1630 gmx::arrayRefFromArray(md->invmass, md->nr),
1639 *ir, md->havePartiallyFrozenAtoms, md->homenr, state, wcycle, constr != nullptr);
1642 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1644 updatePrevStepPullCom(pull_work, state);
1647 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1650 /* ############## IF NOT VV, Calculate globals HERE ############ */
1651 /* With Leap-Frog we can skip compute_globals at
1652 * non-communication steps, but we need to calculate
1653 * the kinetic energy one step before communication.
1656 // Organize to do inter-simulation signalling on steps if
1657 // and when algorithms require it.
1658 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1660 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1662 // Copy coordinates when needed to stop the CM motion.
1663 if (useGpuForUpdate && (bDoReplEx || (!EI_VV(ir->eI) && bStopCM)))
1665 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1666 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1668 // Since we're already communicating at this step, we
1669 // can propagate intra-simulation signals. Note that
1670 // check_nstglobalcomm has the responsibility for
1671 // choosing the value of nstglobalcomm that is one way
1672 // bGStat becomes true, so we can't get into a
1673 // situation where e.g. checkpointing can't be
1675 bool doIntraSimSignal = true;
1676 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1684 makeConstArrayRef(state->x),
1685 makeConstArrayRef(state->v),
1696 (!EI_VV(ir->eI) && bCalcEner && constr != nullptr) ? constr->rmsdData()
1697 : gmx::ArrayRef<real>{},
1701 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1702 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1703 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1704 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1705 | (DOMAINDECOMP(cr) && dd_localTopologyChecker(*cr->dd).shouldCheckNumberOfBondedInteractions()
1706 ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1708 if (DOMAINDECOMP(cr))
1710 dd_localTopologyChecker(cr->dd)->checkNumberOfBondedInteractions(
1711 &top, makeConstArrayRef(state->x), state->box);
1713 if (!EI_VV(ir->eI) && bStopCM)
1715 process_and_stopcm_grp(
1716 fplog, &vcm, *md, makeArrayRef(state->x), makeArrayRef(state->v));
1717 inc_nrnb(nrnb, eNR_STOPCM, md->homenr);
1719 // TODO: The special case of removing CM motion should be dealt more gracefully
1720 if (useGpuForUpdate)
1722 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1723 // Here we block until the H2D copy completes because event sync with the
1724 // force kernels that use the coordinates on the next steps is not implemented
1725 // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1726 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1727 // If the COM removal changed the velocities on the CPU, this has to be accounted for.
1728 if (vcm.mode != ComRemovalAlgorithm::No)
1730 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1737 /* ############# END CALC EKIN AND PRESSURE ################# */
1739 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1740 the virial that should probably be addressed eventually. state->veta has better properies,
1741 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1742 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1744 if (ir->efep != FreeEnergyPerturbationType::No && !EI_VV(ir->eI))
1746 /* Sum up the foreign energy and dK/dl terms for md and sd.
1747 Currently done every step so that dH/dl is correct in the .edr */
1748 accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
1751 update_pcouple_after_coordinates(fplog,
1755 md->cFREEZE ? gmx::arrayRefFromArray(md->cFREEZE, md->nr)
1756 : gmx::ArrayRef<const unsigned short>(),
1766 const bool doBerendsenPressureCoupling = (inputrec->epc == PressureCoupling::Berendsen
1767 && do_per_step(step, inputrec->nstpcouple));
1768 const bool doCRescalePressureCoupling = (inputrec->epc == PressureCoupling::CRescale
1769 && do_per_step(step, inputrec->nstpcouple));
1771 && (doBerendsenPressureCoupling || doCRescalePressureCoupling || doParrinelloRahman))
1773 integrator->scaleCoordinates(pressureCouplingMu);
1774 if (doCRescalePressureCoupling)
1776 matrix pressureCouplingInvMu;
1777 gmx::invertBoxMatrix(pressureCouplingMu, pressureCouplingInvMu);
1778 integrator->scaleVelocities(pressureCouplingInvMu);
1780 integrator->setPbc(PbcType::Xyz, state->box);
1783 /* ################# END UPDATE STEP 2 ################# */
1784 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1786 /* The coordinates (x) were unshifted in update */
1789 /* We will not sum ekinh_old,
1790 * so signal that we still have to do it.
1792 bSumEkinhOld = TRUE;
1797 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1799 /* use the directly determined last velocity, not actually the averaged half steps */
1800 if (bTrotter && ir->eI == IntegrationAlgorithm::VV)
1802 enerd->term[F_EKIN] = last_ekin;
1804 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1806 if (integratorHasConservedEnergyQuantity(ir))
1810 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1814 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1817 /* ######### END PREPARING EDR OUTPUT ########### */
1823 if (fplog && do_log && bDoExpanded)
1825 /* only needed if doing expanded ensemble */
1826 PrintFreeEnergyInfoToFile(fplog,
1828 ir->expandedvals.get(),
1829 ir->bSimTemp ? ir->simtempvals.get() : nullptr,
1830 state_global->dfhist,
1837 energyOutput.addDataAtEnergyStep(bDoDHDL,
1843 ir->expandedvals.get(),
1845 PTCouplingArrays{ state->boxv,
1846 state->nosehoover_xi,
1847 state->nosehoover_vxi,
1849 state->nhpres_vxi },
1859 energyOutput.recordNonEnergyStep();
1862 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1863 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1865 if (doSimulatedAnnealing)
1867 gmx::EnergyOutput::printAnnealingTemperatures(
1868 do_log ? fplog : nullptr, groups, &(ir->opts));
1870 if (do_log || do_ene || do_dr || do_or)
1872 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1876 do_log ? fplog : nullptr,
1882 if (do_log && ir->bDoAwh && awh->hasFepLambdaDimension())
1884 const bool isInitialOutput = false;
1885 printLambdaStateToLog(fplog, state->lambda, isInitialOutput);
1890 pull_print_output(pull_work, step, t);
1893 if (do_per_step(step, ir->nstlog))
1895 if (fflush(fplog) != 0)
1897 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1903 /* Have to do this part _after_ outputting the logfile and the edr file */
1904 /* Gets written into the state at the beginning of next loop*/
1905 state->fep_state = lamnew;
1907 else if (ir->bDoAwh && awh->needForeignEnergyDifferences(step))
1909 state->fep_state = awh->fepLambdaState();
1911 /* Print the remaining wall clock time for the run */
1912 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1916 fprintf(stderr, "\n");
1918 print_time(stderr, walltime_accounting, step, ir, cr);
1921 /* Ion/water position swapping.
1922 * Not done in last step since trajectory writing happens before this call
1923 * in the MD loop and exchanges would be lost anyway. */
1924 bNeedRepartition = FALSE;
1925 if ((ir->eSwapCoords != SwapType::No) && (step > 0) && !bLastStep
1926 && do_per_step(step, ir->swap->nstswap))
1928 bNeedRepartition = do_swapcoords(cr,
1934 as_rvec_array(state->x.data()),
1936 MASTER(cr) && mdrunOptions.verbose,
1939 if (bNeedRepartition && DOMAINDECOMP(cr))
1941 dd_collect_state(cr->dd, state, state_global);
1945 /* Replica exchange */
1949 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1952 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1954 dd_partition_system(fplog,
1975 upd.updateAfterPartition(state->natoms,
1976 md->cFREEZE ? gmx::arrayRefFromArray(md->cFREEZE, md->nr)
1977 : gmx::ArrayRef<const unsigned short>(),
1978 md->cTC ? gmx::arrayRefFromArray(md->cTC, md->nr)
1979 : gmx::ArrayRef<const unsigned short>());
1985 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1986 /* With all integrators, except VV, we need to retain the pressure
1987 * at the current step for coupling at the next step.
1989 if ((state->flags & enumValueToBitMask(StateEntry::PressurePrevious))
1990 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1992 /* Store the pressure in t_state for pressure coupling
1993 * at the next MD step.
1995 copy_mat(pres, state->pres_prev);
1998 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2000 if ((membed != nullptr) && (!bLastStep))
2002 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
2005 cycles = wallcycle_stop(wcycle, WallCycleCounter::Step);
2006 if (DOMAINDECOMP(cr) && wcycle)
2008 dd_cycles_add(cr->dd, cycles, ddCyclStep);
2011 /* increase the MD step number */
2018 fcReportProgress(ir->nsteps + ir->init_step, step);
2022 resetHandler->resetCounters(
2023 step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb, fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
2025 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
2026 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
2028 /* End of main MD loop */
2030 /* Closing TNG files can include compressing data. Therefore it is good to do that
2031 * before stopping the time measurements. */
2032 mdoutf_tng_close(outf);
2034 /* Stop measuring walltime */
2035 walltime_accounting_end_time(walltime_accounting);
2037 if (!thisRankHasDuty(cr, DUTY_PME))
2039 /* Tell the PME only node to finish */
2040 gmx_pme_send_finish(cr);
2045 if (ir->nstcalcenergy > 0)
2047 energyOutput.printEnergyConservation(fplog, ir->simulation_part, EI_MD(ir->eI));
2049 gmx::EnergyOutput::printAnnealingTemperatures(fplog, groups, &(ir->opts));
2050 energyOutput.printAverages(fplog, groups);
2057 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
2060 done_shellfc(fplog, shellfc, step_rel);
2062 if (useReplicaExchange && MASTER(cr))
2064 print_replica_exchange_statistics(fplog, repl_ex);
2067 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
2069 global_stat_destroy(gstat);