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/mdsetup.h"
65 #include "gromacs/domdec/partition.h"
66 #include "gromacs/essentialdynamics/edsam.h"
67 #include "gromacs/ewald/pme_load_balancing.h"
68 #include "gromacs/ewald/pme_pp.h"
69 #include "gromacs/fileio/trxio.h"
70 #include "gromacs/gmxlib/network.h"
71 #include "gromacs/gmxlib/nrnb.h"
72 #include "gromacs/gpu_utils/device_stream_manager.h"
73 #include "gromacs/gpu_utils/gpu_utils.h"
74 #include "gromacs/imd/imd.h"
75 #include "gromacs/listed_forces/listed_forces.h"
76 #include "gromacs/math/functions.h"
77 #include "gromacs/math/invertmatrix.h"
78 #include "gromacs/math/vec.h"
79 #include "gromacs/math/vectypes.h"
80 #include "gromacs/mdlib/checkpointhandler.h"
81 #include "gromacs/mdlib/compute_io.h"
82 #include "gromacs/mdlib/constr.h"
83 #include "gromacs/mdlib/coupling.h"
84 #include "gromacs/mdlib/ebin.h"
85 #include "gromacs/mdlib/enerdata_utils.h"
86 #include "gromacs/mdlib/energyoutput.h"
87 #include "gromacs/mdlib/expanded.h"
88 #include "gromacs/mdlib/force.h"
89 #include "gromacs/mdlib/force_flags.h"
90 #include "gromacs/mdlib/forcerec.h"
91 #include "gromacs/mdlib/freeenergyparameters.h"
92 #include "gromacs/mdlib/md_support.h"
93 #include "gromacs/mdlib/mdatoms.h"
94 #include "gromacs/mdlib/mdoutf.h"
95 #include "gromacs/mdlib/membed.h"
96 #include "gromacs/mdlib/resethandler.h"
97 #include "gromacs/mdlib/sighandler.h"
98 #include "gromacs/mdlib/simulationsignal.h"
99 #include "gromacs/mdlib/stat.h"
100 #include "gromacs/mdlib/stophandler.h"
101 #include "gromacs/mdlib/tgroup.h"
102 #include "gromacs/mdlib/trajectory_writing.h"
103 #include "gromacs/mdlib/update.h"
104 #include "gromacs/mdlib/update_constrain_gpu.h"
105 #include "gromacs/mdlib/update_vv.h"
106 #include "gromacs/mdlib/vcm.h"
107 #include "gromacs/mdlib/vsite.h"
108 #include "gromacs/mdrunutility/handlerestart.h"
109 #include "gromacs/mdrunutility/multisim.h"
110 #include "gromacs/mdrunutility/printtime.h"
111 #include "gromacs/mdtypes/awh_history.h"
112 #include "gromacs/mdtypes/awh_params.h"
113 #include "gromacs/mdtypes/commrec.h"
114 #include "gromacs/mdtypes/df_history.h"
115 #include "gromacs/mdtypes/energyhistory.h"
116 #include "gromacs/mdtypes/fcdata.h"
117 #include "gromacs/mdtypes/forcebuffers.h"
118 #include "gromacs/mdtypes/forcerec.h"
119 #include "gromacs/mdtypes/group.h"
120 #include "gromacs/mdtypes/inputrec.h"
121 #include "gromacs/mdtypes/interaction_const.h"
122 #include "gromacs/mdtypes/md_enums.h"
123 #include "gromacs/mdtypes/mdatom.h"
124 #include "gromacs/mdtypes/mdrunoptions.h"
125 #include "gromacs/mdtypes/multipletimestepping.h"
126 #include "gromacs/mdtypes/observableshistory.h"
127 #include "gromacs/mdtypes/pullhistory.h"
128 #include "gromacs/mdtypes/simulation_workload.h"
129 #include "gromacs/mdtypes/state.h"
130 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
131 #include "gromacs/modularsimulator/energydata.h"
132 #include "gromacs/nbnxm/gpu_data_mgmt.h"
133 #include "gromacs/nbnxm/nbnxm.h"
134 #include "gromacs/pbcutil/pbc.h"
135 #include "gromacs/pulling/output.h"
136 #include "gromacs/pulling/pull.h"
137 #include "gromacs/swap/swapcoords.h"
138 #include "gromacs/timing/wallcycle.h"
139 #include "gromacs/timing/walltime_accounting.h"
140 #include "gromacs/topology/atoms.h"
141 #include "gromacs/topology/idef.h"
142 #include "gromacs/topology/mtop_util.h"
143 #include "gromacs/topology/topology.h"
144 #include "gromacs/trajectory/trajectoryframe.h"
145 #include "gromacs/utility/basedefinitions.h"
146 #include "gromacs/utility/cstringutil.h"
147 #include "gromacs/utility/fatalerror.h"
148 #include "gromacs/utility/logger.h"
149 #include "gromacs/utility/real.h"
150 #include "gromacs/utility/smalloc.h"
152 #include "legacysimulator.h"
153 #include "replicaexchange.h"
156 using gmx::SimulationSignaller;
158 void gmx::LegacySimulator::do_md()
160 // TODO Historically, the EM and MD "integrators" used different
161 // names for the t_inputrec *parameter, but these must have the
162 // same name, now that it's a member of a struct. We use this ir
163 // alias to avoid a large ripple of nearly useless changes.
164 // t_inputrec is being replaced by IMdpOptionsProvider, so this
165 // will go away eventually.
166 const t_inputrec* ir = inputrec;
168 int64_t step, step_rel;
169 double t, t0 = ir->init_t;
170 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
171 gmx_bool bNS = FALSE, bNStList, bStopCM, bFirstStep, bInitStep, bLastStep = FALSE;
172 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
173 gmx_bool do_ene, do_log, do_verbose;
174 gmx_bool bMasterState;
175 unsigned int force_flags;
176 tensor force_vir = { { 0 } }, shake_vir = { { 0 } }, total_vir = { { 0 } }, pres = { { 0 } };
179 matrix pressureCouplingMu, M;
180 gmx_repl_ex_t repl_ex = nullptr;
181 gmx_global_stat_t gstat;
182 gmx_shellfc_t* shellfc;
183 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
186 std::vector<RVec> cbuf;
192 real saved_conserved_quantity = 0;
195 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
197 /* PME load balancing data for GPU kernels */
198 gmx_bool bPMETune = FALSE;
199 gmx_bool bPMETunePrinting = FALSE;
201 bool bInteractiveMDstep = false;
203 /* Domain decomposition could incorrectly miss a bonded
204 interaction, but checking for that requires a global
205 communication stage, which does not otherwise happen in DD
206 code. So we do that alongside the first global energy reduction
207 after a new DD is made. These variables handle whether the
208 check happens, and the result it returns. */
209 bool shouldCheckNumberOfBondedInteractions = false;
210 int totalNumberOfBondedInteractions = -1;
212 SimulationSignals signals;
213 // Most global communnication stages don't propagate mdrun
214 // signals, and will use this object to achieve that.
215 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
217 if (!mdrunOptions.writeConfout)
219 // This is on by default, and the main known use case for
220 // turning it off is for convenience in benchmarking, which is
221 // something that should not show up in the general user
226 "The -noconfout functionality is deprecated, and may be removed in a "
230 /* md-vv uses averaged full step velocities for T-control
231 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
232 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
233 bTrotter = (EI_VV(ir->eI)
234 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
236 const bool bRerunMD = false;
238 int nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
239 bGStatEveryStep = (nstglobalcomm == 1);
241 const SimulationGroups* groups = &top_global.groups;
243 std::unique_ptr<EssentialDynamics> ed = nullptr;
244 if (opt2bSet("-ei", nfile, fnm))
246 /* Initialize essential dynamics sampling */
247 ed = init_edsam(mdlog,
248 opt2fn_null("-ei", nfile, fnm),
249 opt2fn("-eo", nfile, fnm),
259 else if (observablesHistory->edsamHistory)
262 "The checkpoint is from a run with essential dynamics sampling, "
263 "but the current run did not specify the -ei option. "
264 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
267 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
268 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
269 initialize_lambdas(fplog, *ir, MASTER(cr), fep_state, lambda);
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 =
288 (fcdata.disres->nsystems > 1) || ((ms != nullptr) && (fcdata.orires->nr != 0));
289 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim() && (ms != nullptr));
291 // Replica exchange, ensemble restraints and AWH need all
292 // simulations to remain synchronized, so they need
293 // checkpoints and stop conditions to act on the same step, so
294 // the propagation of such signals must take place between
295 // simulations, not just within simulations.
296 // TODO: Make algorithm initializers set these flags.
297 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
299 if (simulationsShareState)
301 // Inter-simulation signal communication does not need to happen
302 // often, so we use a minimum of 200 steps to reduce overhead.
303 const int c_minimumInterSimulationSignallingInterval = 200;
304 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
309 if (startingBehavior != StartingBehavior::RestartWithAppending)
311 pleaseCiteCouplingAlgorithms(fplog, *ir);
313 gmx_mdoutf* outf = init_mdoutf(fplog,
325 simulationsShareState,
327 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
331 mdoutf_get_fp_dhdl(outf),
334 simulationsShareState,
337 gstat = global_stat_init(ir);
339 const auto& simulationWork = runScheduleWork->simulationWork;
340 const bool useGpuForPme = simulationWork.useGpuPme;
341 const bool useGpuForNonbonded = simulationWork.useGpuNonbonded;
342 const bool useGpuForBufferOps = simulationWork.useGpuBufferOps;
343 const bool useGpuForUpdate = simulationWork.useGpuUpdate;
345 /* Check for polarizable models and flexible constraints */
346 shellfc = init_shell_flexcon(fplog,
348 constr ? constr->numFlexibleConstraints() : 0,
354 double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
355 if ((io > 2000) && MASTER(cr))
357 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
361 // Local state only becomes valid now.
362 std::unique_ptr<t_state> stateInstance;
365 gmx_localtop_t top(top_global.ffparams);
367 auto mdatoms = mdAtoms->mdatoms();
369 ForceBuffers f(fr->useMts,
370 ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
371 ? PinningPolicy::PinnedIfSupported
372 : PinningPolicy::CannotBePinned);
373 if (DOMAINDECOMP(cr))
375 stateInstance = std::make_unique<t_state>();
376 state = stateInstance.get();
377 dd_init_local_state(*cr->dd, state_global, state);
379 /* Distribute the charge groups over the nodes from the master node */
380 dd_partition_system(fplog,
401 shouldCheckNumberOfBondedInteractions = true;
402 upd.setNumAtoms(state->natoms);
406 state_change_natoms(state_global, state_global->natoms);
407 /* Copy the pointer to the global state */
408 state = state_global;
410 /* Generate and initialize new topology */
411 mdAlgorithmsSetupAtomData(cr, *ir, top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
413 upd.setNumAtoms(state->natoms);
416 std::unique_ptr<UpdateConstrainGpu> integrator;
418 StatePropagatorDataGpu* stateGpu = fr->stateGpu;
420 // TODO: the assertions below should be handled by UpdateConstraintsBuilder.
423 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
424 || constr->numConstraintsTotal() == 0,
425 "Constraints in domain decomposition are only supported with update "
426 "groups if using GPU update.\n");
427 GMX_RELEASE_ASSERT(ir->eConstrAlg != ConstraintAlgorithm::Shake || constr == nullptr
428 || constr->numConstraintsTotal() == 0,
429 "SHAKE is not supported with GPU update.");
430 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
431 "Either PME or short-ranged non-bonded interaction tasks must run on "
432 "the GPU to use GPU update.\n");
433 GMX_RELEASE_ASSERT(ir->eI == IntegrationAlgorithm::MD,
434 "Only the md integrator is supported with the GPU update.\n");
436 ir->etc != TemperatureCoupling::NoseHoover,
437 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
439 ir->epc == PressureCoupling::No || ir->epc == PressureCoupling::ParrinelloRahman
440 || ir->epc == PressureCoupling::Berendsen || ir->epc == PressureCoupling::CRescale,
441 "Only Parrinello-Rahman, Berendsen, and C-rescale pressure coupling are supported "
442 "with the GPU update.\n");
443 GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
444 "Virtual sites are not supported with the GPU update.\n");
445 GMX_RELEASE_ASSERT(ed == nullptr,
446 "Essential dynamics is not supported with the GPU update.\n");
447 GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(*ir->pull),
448 "Constraints pulling is not supported with the GPU update.\n");
449 GMX_RELEASE_ASSERT(fcdata.orires->nr == 0,
450 "Orientation restraints are not supported with the GPU update.\n");
452 ir->efep == FreeEnergyPerturbationType::No
453 || (!haveFepPerturbedMasses(top_global) && !havePerturbedConstraints(top_global)),
454 "Free energy perturbation of masses and constraints are not supported with the GPU "
457 if (constr != nullptr && constr->numConstraintsTotal() > 0)
461 .appendText("Updating coordinates and applying constraints on the GPU.");
465 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
467 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
468 "Device stream manager should be initialized in order to use GPU "
469 "update-constraints.");
471 fr->deviceStreamManager->streamIsValid(gmx::DeviceStreamType::UpdateAndConstraints),
472 "Update stream should be initialized in order to use GPU "
473 "update-constraints.");
474 integrator = std::make_unique<UpdateConstrainGpu>(
478 fr->deviceStreamManager->context(),
479 fr->deviceStreamManager->stream(gmx::DeviceStreamType::UpdateAndConstraints),
480 stateGpu->xUpdatedOnDevice(),
483 integrator->setPbc(PbcType::Xyz, state->box);
486 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
488 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
492 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
495 // NOTE: The global state is no longer used at this point.
496 // But state_global is still used as temporary storage space for writing
497 // the global state to file and potentially for replica exchange.
498 // (Global topology should persist.)
500 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
504 /* Check nstexpanded here, because the grompp check was broken */
505 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
508 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
510 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
515 EnergyData::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
518 preparePrevStepPullCom(ir,
520 gmx::arrayRefFromArray(mdatoms->massT, mdatoms->nr),
524 startingBehavior != StartingBehavior::NewSimulation);
526 // TODO: Remove this by converting AWH into a ForceProvider
527 auto awh = prepareAwhModule(fplog,
532 startingBehavior != StartingBehavior::NewSimulation,
534 opt2fn("-awh", nfile, fnm),
537 if (useReplicaExchange && MASTER(cr))
539 repl_ex = init_replica_exchange(fplog, ms, top_global.natoms, ir, replExParams);
541 /* PME tuning is only supported in the Verlet scheme, with PME for
542 * Coulomb. It is not supported with only LJ PME. */
543 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
544 && ir->cutoff_scheme != CutoffScheme::Group);
546 pme_load_balancing_t* pme_loadbal = nullptr;
550 &pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata, fr->nbv->useGpu());
553 if (!ir->bContinuation)
555 if (state->flags & enumValueToBitMask(StateEntry::V))
557 auto v = makeArrayRef(state->v);
558 /* Set the velocities of vsites, shells and frozen atoms to zero */
559 for (i = 0; i < mdatoms->homenr; i++)
561 if (mdatoms->ptype[i] == ParticleType::Shell)
565 else if (mdatoms->cFREEZE)
567 for (m = 0; m < DIM; m++)
569 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
580 /* Constrain the initial coordinates and velocities */
581 do_constrain_first(fplog,
586 state->x.arrayRefWithPadding(),
587 state->v.arrayRefWithPadding(),
589 state->lambda[FreeEnergyPerturbationCouplingType::Bonded]);
593 if (ir->efep != FreeEnergyPerturbationType::No)
595 /* Set free energy calculation frequency as the greatest common
596 * denominator of nstdhdl and repl_ex_nst. */
597 nstfep = ir->fepvals->nstdhdl;
600 nstfep = std::gcd(ir->expandedvals->nstexpanded, nstfep);
602 if (useReplicaExchange)
604 nstfep = std::gcd(replExParams.exchangeInterval, nstfep);
608 nstfep = std::gcd(ir->awhParams->nstSampleCoord(), nstfep);
612 /* Be REALLY careful about what flags you set here. You CANNOT assume
613 * this is the first step, since we might be restarting from a checkpoint,
614 * and in that case we should not do any modifications to the state.
616 bStopCM = (ir->comm_mode != ComRemovalAlgorithm::No && !ir->bContinuation);
618 // When restarting from a checkpoint, it can be appropriate to
619 // initialize ekind from quantities in the checkpoint. Otherwise,
620 // compute_globals must initialize ekind before the simulation
621 // starts/restarts. However, only the master rank knows what was
622 // found in the checkpoint file, so we have to communicate in
623 // order to coordinate the restart.
625 // TODO Consider removing this communication if/when checkpoint
626 // reading directly follows .tpr reading, because all ranks can
627 // agree on hasReadEkinState at that time.
628 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
631 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr->mpi_comm_mygroup);
633 if (hasReadEkinState)
635 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
638 unsigned int cglo_flags =
639 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
640 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
642 bSumEkinhOld = FALSE;
644 t_vcm vcm(top_global.groups, *ir);
645 reportComRemovalInfo(fplog, vcm);
647 /* To minimize communication, compute_globals computes the COM velocity
648 * and the kinetic energy for the velocities without COM motion removed.
649 * Thus to get the kinetic energy without the COM contribution, we need
650 * to call compute_globals twice.
652 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
654 unsigned int cglo_flags_iteration = cglo_flags;
655 if (bStopCM && cgloIteration == 0)
657 cglo_flags_iteration |= CGLO_STOPCM;
658 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
660 compute_globals(gstat,
665 makeConstArrayRef(state->x),
666 makeConstArrayRef(state->v),
677 gmx::ArrayRef<real>{},
680 &totalNumberOfBondedInteractions,
683 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
685 if (cglo_flags_iteration & CGLO_STOPCM)
687 /* At initialization, do not pass x with acceleration-correction mode
688 * to avoid (incorrect) correction of the initial coordinates.
690 auto x = (vcm.mode == ComRemovalAlgorithm::LinearAccelerationCorrection)
692 : makeArrayRef(state->x);
693 process_and_stopcm_grp(fplog, &vcm, *mdatoms, x, makeArrayRef(state->v));
694 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
697 checkNumberOfBondedInteractions(mdlog,
699 totalNumberOfBondedInteractions,
702 makeConstArrayRef(state->x),
704 &shouldCheckNumberOfBondedInteractions);
705 if (ir->eI == IntegrationAlgorithm::VVAK)
707 /* a second call to get the half step temperature initialized as well */
708 /* we do the same call as above, but turn the pressure off -- internally to
709 compute_globals, this is recognized as a velocity verlet half-step
710 kinetic energy calculation. This minimized excess variables, but
711 perhaps loses some logic?*/
713 compute_globals(gstat,
718 makeConstArrayRef(state->x),
719 makeConstArrayRef(state->v),
730 gmx::ArrayRef<real>{},
735 cglo_flags & ~CGLO_PRESSURE);
738 /* Calculate the initial half step temperature, and save the ekinh_old */
739 if (startingBehavior == StartingBehavior::NewSimulation)
741 for (i = 0; (i < ir->opts.ngtc); i++)
743 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
747 /* need to make an initiation call to get the Trotter variables set, as well as other constants
748 for non-trotter temperature control */
749 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
753 if (!ir->bContinuation)
755 if (constr && ir->eConstrAlg == ConstraintAlgorithm::Lincs)
758 "RMS relative constraint deviation after constraining: %.2e\n",
761 if (EI_STATE_VELOCITY(ir->eI))
763 real temp = enerd->term[F_TEMP];
764 if (ir->eI != IntegrationAlgorithm::VV)
766 /* Result of Ekin averaged over velocities of -half
767 * and +half step, while we only have -half step here.
771 fprintf(fplog, "Initial temperature: %g K\n", temp);
776 fprintf(stderr, "starting mdrun '%s'\n", *(top_global.name));
779 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
783 sprintf(tbuf, "%s", "infinite");
785 if (ir->init_step > 0)
788 "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
789 gmx_step_str(ir->init_step + ir->nsteps, sbuf),
791 gmx_step_str(ir->init_step, sbuf2),
792 ir->init_step * ir->delta_t);
796 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
798 fprintf(fplog, "\n");
801 walltime_accounting_start_time(walltime_accounting);
802 wallcycle_start(wcycle, ewcRUN);
803 print_start(fplog, cr, walltime_accounting, "mdrun");
805 /***********************************************************
809 ************************************************************/
812 /* Skip the first Nose-Hoover integration when we get the state from tpx */
813 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
814 bSumEkinhOld = FALSE;
816 bNeedRepartition = FALSE;
818 step = ir->init_step;
821 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
822 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
823 simulationsShareState,
826 mdrunOptions.reproducible,
828 mdrunOptions.maximumHoursToRun,
833 walltime_accounting);
835 auto checkpointHandler = std::make_unique<CheckpointHandler>(
836 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]),
837 simulationsShareState,
840 mdrunOptions.writeConfout,
841 mdrunOptions.checkpointOptions.period);
843 const bool resetCountersIsLocal = true;
844 auto resetHandler = std::make_unique<ResetHandler>(
845 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
846 !resetCountersIsLocal,
849 mdrunOptions.timingOptions.resetHalfway,
850 mdrunOptions.maximumHoursToRun,
853 walltime_accounting);
855 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
857 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
859 logInitialMultisimStatus(ms, cr, mdlog, simulationsShareState, ir->nsteps, ir->init_step);
862 /* and stop now if we should */
863 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
867 /* Determine if this is a neighbor search step */
868 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
870 if (bPMETune && bNStList)
872 // This has to be here because PME load balancing is called so early.
873 // TODO: Move to after all booleans are defined.
874 if (useGpuForUpdate && !bFirstStep)
876 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
877 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
879 /* PME grid + cut-off optimization with GPUs or PME nodes */
880 pme_loadbal_do(pme_loadbal,
882 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
893 simulationWork.useGpuPmePpCommunication);
896 wallcycle_start(wcycle, ewcSTEP);
898 bLastStep = (step_rel == ir->nsteps);
899 t = t0 + step * ir->delta_t;
901 // TODO Refactor this, so that nstfep does not need a default value of zero
902 if (ir->efep != FreeEnergyPerturbationType::No || ir->bSimTemp)
904 /* find and set the current lambdas */
905 state->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
907 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
908 bDoFEP = ((ir->efep != FreeEnergyPerturbationType::No) && do_per_step(step, nstfep));
909 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
913 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
914 && do_per_step(step, replExParams.exchangeInterval));
916 if (doSimulatedAnnealing)
918 // TODO: Avoid changing inputrec (#3854)
919 // Simulated annealing updates the reference temperature.
920 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
921 update_annealing_target_temp(nonConstInputrec, t, &upd);
924 /* Stop Center of Mass motion */
925 bStopCM = (ir->comm_mode != ComRemovalAlgorithm::No && do_per_step(step, ir->nstcomm));
927 /* Determine whether or not to do Neighbour Searching */
928 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
930 /* Note that the stopHandler will cause termination at nstglobalcomm
931 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
932 * nstpcouple steps, we have computed the half-step kinetic energy
933 * of the previous step and can always output energies at the last step.
935 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
937 /* do_log triggers energy and virial calculation. Because this leads
938 * to different code paths, forces can be different. Thus for exact
939 * continuation we should avoid extra log output.
940 * Note that the || bLastStep can result in non-exact continuation
941 * beyond the last step. But we don't consider that to be an issue.
943 do_log = (do_per_step(step, ir->nstlog)
944 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
945 do_verbose = mdrunOptions.verbose
946 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
948 if (useGpuForUpdate && !bFirstStep && bNS)
950 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
951 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
952 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
953 // Copy coordinate from the GPU when needed at the search step.
954 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
955 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
956 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
957 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
960 // We only need to calculate virtual velocities if we are writing them in the current step
961 const bool needVirtualVelocitiesThisStep =
963 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep());
965 if (vsite != nullptr)
967 // Virtual sites need to be updated before domain decomposition and forces are calculated
968 wallcycle_start(wcycle, ewcVSITECONSTR);
969 // md-vv calculates virtual velocities once it has full-step real velocities
970 vsite->construct(state->x,
973 (!EI_VV(inputrec->eI) && needVirtualVelocitiesThisStep)
974 ? VSiteOperation::PositionsAndVelocities
975 : VSiteOperation::Positions);
976 wallcycle_stop(wcycle, ewcVSITECONSTR);
979 if (bNS && !(bFirstStep && ir->bContinuation))
981 bMasterState = FALSE;
982 /* Correct the new box if it is too skewed */
983 if (inputrecDynamicBox(ir))
985 if (correct_box(fplog, step, state->box))
988 // If update is offloaded, it should be informed about the box size change
991 integrator->setPbc(PbcType::Xyz, state->box);
995 if (DOMAINDECOMP(cr) && bMasterState)
997 dd_collect_state(cr->dd, state, state_global);
1000 if (DOMAINDECOMP(cr))
1002 /* Repartition the domain decomposition */
1003 dd_partition_system(fplog,
1023 do_verbose && !bPMETunePrinting);
1024 shouldCheckNumberOfBondedInteractions = true;
1025 upd.setNumAtoms(state->natoms);
1029 // Allocate or re-size GPU halo exchange object, if necessary
1030 if (bNS && havePPDomainDecomposition(cr) && simulationWork.useGpuHaloExchange)
1032 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
1033 "GPU device manager has to be initialized to use GPU "
1034 "version of halo exchange.");
1035 constructGpuHaloExchange(mdlog, *cr, *fr->deviceStreamManager, wcycle);
1038 if (MASTER(cr) && do_log)
1040 gmx::EnergyOutput::printHeader(
1041 fplog, step, t); /* can we improve the information printed here? */
1044 if (ir->efep != FreeEnergyPerturbationType::No)
1046 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
1051 /* We need the kinetic energy at minus the half step for determining
1052 * the full step kinetic energy and possibly for T-coupling.*/
1053 /* This may not be quite working correctly yet . . . . */
1054 compute_globals(gstat,
1059 makeConstArrayRef(state->x),
1060 makeConstArrayRef(state->v),
1071 gmx::ArrayRef<real>{},
1074 &totalNumberOfBondedInteractions,
1076 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1077 checkNumberOfBondedInteractions(mdlog,
1079 totalNumberOfBondedInteractions,
1082 makeConstArrayRef(state->x),
1084 &shouldCheckNumberOfBondedInteractions);
1086 clear_mat(force_vir);
1088 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
1090 /* Determine the energy and pressure:
1091 * at nstcalcenergy steps and at energy output steps (set below).
1093 if (EI_VV(ir->eI) && (!bInitStep))
1095 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1096 bCalcVir = bCalcEnerStep
1097 || (ir->epc != PressureCoupling::No
1098 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
1102 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1103 bCalcVir = bCalcEnerStep
1104 || (ir->epc != PressureCoupling::No && do_per_step(step, ir->nstpcouple));
1106 bCalcEner = bCalcEnerStep;
1108 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1110 if (do_ene || do_log || bDoReplEx)
1116 /* Do we need global communication ? */
1117 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
1118 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
1120 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
1121 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
1122 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
1123 if (fr->useMts && !do_per_step(step, ir->nstfout))
1125 force_flags |= GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE;
1130 /* Now is the time to relax the shells */
1131 relax_shell_flexcon(fplog,
1134 mdrunOptions.verbose,
1146 state->x.arrayRefWithPadding(),
1147 state->v.arrayRefWithPadding(),
1162 ddBalanceRegionHandler);
1166 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
1167 is updated (or the AWH update will be performed twice for one step when continuing).
1168 It would be best to call this update function from do_md_trajectory_writing but that
1169 would occur after do_force. One would have to divide the update_awh function into one
1170 function applying the AWH force and one doing the AWH bias update. The update AWH
1171 bias function could then be called after do_md_trajectory_writing (then containing
1172 update_awh_history). The checkpointing will in the future probably moved to the start
1173 of the md loop which will rid of this issue. */
1174 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
1176 awh->updateHistory(state_global->awhHistory.get());
1179 /* The coordinates (x) are shifted (to get whole molecules)
1181 * This is parallellized as well, and does communication too.
1182 * Check comments in sim_util.c
1197 state->x.arrayRefWithPadding(),
1209 ed ? ed->getLegacyED() : nullptr,
1210 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1211 ddBalanceRegionHandler);
1214 // VV integrators do not need the following velocity half step
1215 // if it is the first step after starting from a checkpoint.
1216 // That is, the half step is needed on all other steps, and
1217 // also the first step when starting from a .tpr file.
1220 integrateVVFirstStep(step,
1253 &shouldCheckNumberOfBondedInteractions,
1254 &saved_conserved_quantity,
1264 if (vsite != nullptr && needVirtualVelocitiesThisStep)
1266 // Positions were calculated earlier
1267 wallcycle_start(wcycle, ewcVSITECONSTR);
1268 vsite->construct(state->x, state->v, state->box, VSiteOperation::Velocities);
1269 wallcycle_stop(wcycle, ewcVSITECONSTR);
1273 /* ######## END FIRST UPDATE STEP ############## */
1274 /* ######## If doing VV, we now have v(dt) ###### */
1277 /* perform extended ensemble sampling in lambda - we don't
1278 actually move to the new state before outputting
1279 statistics, but if performing simulated tempering, we
1280 do update the velocities and the tau_t. */
1281 // TODO: Avoid changing inputrec (#3854)
1282 // Simulated tempering updates the reference temperature.
1283 // Expanded ensemble without simulated tempering does not change the inputrec.
1284 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
1285 lamnew = ExpandedEnsembleDynamics(fplog,
1293 state->v.rvec_array(),
1295 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1298 copy_df_history(state_global->dfhist, state->dfhist);
1302 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1303 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1304 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1305 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1306 || checkpointHandler->isCheckpointingStep()))
1308 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1309 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1311 // Copy velocities if needed for the output/checkpointing.
1312 // NOTE: Copy on the search steps is done at the beginning of the step.
1313 if (useGpuForUpdate && !bNS
1314 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1316 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1317 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1319 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1320 // and update is offloaded hence forces are kept on the GPU for update and have not been
1321 // already transferred in do_force().
1322 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1323 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1324 // prior to GPU update.
1325 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1326 // copy call in do_force(...).
1327 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1328 // on host after the D2H copy in do_force(...).
1329 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1330 && do_per_step(step, ir->nstfout))
1332 stateGpu->copyForcesFromGpu(f.view().force(), AtomLocality::Local);
1333 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1335 /* Now we have the energies and forces corresponding to the
1336 * coordinates at time t. We must output all of this before
1339 do_md_trajectory_writing(fplog,
1356 checkpointHandler->isCheckpointingStep(),
1359 mdrunOptions.writeConfout,
1361 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1362 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1364 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1365 if (startingBehavior != StartingBehavior::NewSimulation && bFirstStep
1366 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1368 copy_mat(state->svir_prev, shake_vir);
1369 copy_mat(state->fvir_prev, force_vir);
1372 stopHandler->setSignal();
1373 resetHandler->setSignal(walltime_accounting);
1375 if (bGStat || !PAR(cr))
1377 /* In parallel we only have to check for checkpointing in steps
1378 * where we do global communication,
1379 * otherwise the other nodes don't know.
1381 checkpointHandler->setSignal(walltime_accounting);
1384 /* ######### START SECOND UPDATE STEP ################# */
1386 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1387 controlled in preprocessing */
1389 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1391 gmx_bool bIfRandomize;
1392 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1393 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1394 if (constr && bIfRandomize)
1396 constrain_velocities(constr, do_log, do_ene, step, state, nullptr, false, nullptr);
1399 /* Box is changed in update() when we do pressure coupling,
1400 * but we should still use the old box for energy corrections and when
1401 * writing it to the energy file, so it matches the trajectory files for
1402 * the same timestep above. Make a copy in a separate array.
1404 copy_mat(state->box, lastbox);
1408 if (!useGpuForUpdate)
1410 wallcycle_start(wcycle, ewcUPDATE);
1412 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1415 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1416 /* We can only do Berendsen coupling after we have summed
1417 * the kinetic energy or virial. Since the happens
1418 * in global_state after update, we should only do it at
1419 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1424 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1425 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1428 /* With leap-frog type integrators we compute the kinetic energy
1429 * at a whole time step as the average of the half-time step kinetic
1430 * energies of two subsequent steps. Therefore we need to compute the
1431 * half step kinetic energy also if we need energies at the next step.
1433 const bool needHalfStepKineticEnergy =
1434 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1436 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1437 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1438 const bool doParrinelloRahman = (ir->epc == PressureCoupling::ParrinelloRahman
1439 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1443 GMX_ASSERT(!useGpuForUpdate, "GPU update is not supported with VVAK integrator.");
1445 integrateVVSecondStep(step,
1481 if (useGpuForUpdate)
1483 if (bNS && (bFirstStep || DOMAINDECOMP(cr)))
1485 integrator->set(stateGpu->getCoordinates(),
1486 stateGpu->getVelocities(),
1487 stateGpu->getForces(),
1491 // Copy data to the GPU after buffers might have being reinitialized
1492 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1493 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1496 if (simulationWork.useGpuPme && !runScheduleWork->simulationWork.useGpuPmePpCommunication
1497 && !thisRankHasDuty(cr, DUTY_PME))
1499 // The PME forces were recieved to the host, so have to be copied
1500 stateGpu->copyForcesToGpu(f.view().force(), AtomLocality::All);
1502 else if (!runScheduleWork->stepWork.useGpuFBufferOps)
1504 // The buffer ops were not offloaded this step, so the forces are on the
1505 // host and have to be copied
1506 stateGpu->copyForcesToGpu(f.view().force(), AtomLocality::Local);
1509 const bool doTemperatureScaling =
1510 (ir->etc != TemperatureCoupling::No
1511 && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1513 // This applies Leap-Frog, LINCS and SETTLE in succession
1514 integrator->integrate(
1515 stateGpu->getForcesReadyOnDeviceEvent(
1516 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1521 doTemperatureScaling,
1524 ir->nstpcouple * ir->delta_t,
1527 // Copy velocities D2H after update if:
1528 // - Globals are computed this step (includes the energy output steps).
1529 // - Temperature is needed for the next step.
1530 if (bGStat || needHalfStepKineticEnergy)
1532 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1533 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1538 /* With multiple time stepping we need to do an additional normal
1539 * update step to obtain the virial, as the actual MTS integration
1540 * using an acceleration where the slow forces are multiplied by mtsFactor.
1541 * Using that acceleration would result in a virial with the slow
1542 * force contribution would be a factor mtsFactor too large.
1544 if (fr->useMts && bCalcVir && constr != nullptr)
1546 upd.update_for_constraint_virial(
1547 *ir, *mdatoms, *state, f.view().forceWithPadding(), *ekind);
1549 constrain_coordinates(constr,
1554 upd.xp()->arrayRefWithPadding(),
1560 ArrayRefWithPadding<const RVec> forceCombined =
1561 (fr->useMts && step % ir->mtsLevels[1].stepFactor == 0)
1562 ? f.view().forceMtsCombinedWithPadding()
1563 : f.view().forceWithPadding();
1565 *ir, step, mdatoms, state, forceCombined, fcdata, ekind, M, etrtPOSITION, cr, constr != nullptr);
1567 wallcycle_stop(wcycle, ewcUPDATE);
1569 constrain_coordinates(constr,
1574 upd.xp()->arrayRefWithPadding(),
1576 bCalcVir && !fr->useMts,
1579 upd.update_sd_second_half(
1580 *ir, step, &dvdl_constr, mdatoms, state, cr, nrnb, wcycle, constr, do_log, do_ene);
1581 upd.finish_update(*ir, mdatoms, state, wcycle, constr != nullptr);
1584 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1586 updatePrevStepPullCom(pull_work, state);
1589 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1592 /* ############## IF NOT VV, Calculate globals HERE ############ */
1593 /* With Leap-Frog we can skip compute_globals at
1594 * non-communication steps, but we need to calculate
1595 * the kinetic energy one step before communication.
1598 // Organize to do inter-simulation signalling on steps if
1599 // and when algorithms require it.
1600 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1602 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1604 // Copy coordinates when needed to stop the CM motion.
1605 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1607 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1608 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1610 // Since we're already communicating at this step, we
1611 // can propagate intra-simulation signals. Note that
1612 // check_nstglobalcomm has the responsibility for
1613 // choosing the value of nstglobalcomm that is one way
1614 // bGStat becomes true, so we can't get into a
1615 // situation where e.g. checkpointing can't be
1617 bool doIntraSimSignal = true;
1618 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1626 makeConstArrayRef(state->x),
1627 makeConstArrayRef(state->v),
1638 (!EI_VV(ir->eI) && bCalcEner && constr != nullptr) ? constr->rmsdData()
1639 : gmx::ArrayRef<real>{},
1642 &totalNumberOfBondedInteractions,
1644 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1645 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1646 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1647 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1648 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1650 checkNumberOfBondedInteractions(mdlog,
1652 totalNumberOfBondedInteractions,
1655 makeConstArrayRef(state->x),
1657 &shouldCheckNumberOfBondedInteractions);
1658 if (!EI_VV(ir->eI) && bStopCM)
1660 process_and_stopcm_grp(
1661 fplog, &vcm, *mdatoms, makeArrayRef(state->x), makeArrayRef(state->v));
1662 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1664 // TODO: The special case of removing CM motion should be dealt more gracefully
1665 if (useGpuForUpdate)
1667 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1668 // Here we block until the H2D copy completes because event sync with the
1669 // force kernels that use the coordinates on the next steps is not implemented
1670 // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1671 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1672 // If the COM removal changed the velocities on the CPU, this has to be accounted for.
1673 if (vcm.mode != ComRemovalAlgorithm::No)
1675 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1682 /* ############# END CALC EKIN AND PRESSURE ################# */
1684 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1685 the virial that should probably be addressed eventually. state->veta has better properies,
1686 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1687 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1689 if (ir->efep != FreeEnergyPerturbationType::No && !EI_VV(ir->eI))
1691 /* Sum up the foreign energy and dK/dl terms for md and sd.
1692 Currently done every step so that dH/dl is correct in the .edr */
1693 accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
1696 update_pcouple_after_coordinates(
1697 fplog, step, ir, mdatoms, pres, force_vir, shake_vir, pressureCouplingMu, state, nrnb, upd.deform(), !useGpuForUpdate);
1699 const bool doBerendsenPressureCoupling = (inputrec->epc == PressureCoupling::Berendsen
1700 && do_per_step(step, inputrec->nstpcouple));
1701 const bool doCRescalePressureCoupling = (inputrec->epc == PressureCoupling::CRescale
1702 && do_per_step(step, inputrec->nstpcouple));
1704 && (doBerendsenPressureCoupling || doCRescalePressureCoupling || doParrinelloRahman))
1706 integrator->scaleCoordinates(pressureCouplingMu);
1707 if (doCRescalePressureCoupling)
1709 matrix pressureCouplingInvMu;
1710 gmx::invertBoxMatrix(pressureCouplingMu, pressureCouplingInvMu);
1711 integrator->scaleVelocities(pressureCouplingInvMu);
1713 integrator->setPbc(PbcType::Xyz, state->box);
1716 /* ################# END UPDATE STEP 2 ################# */
1717 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1719 /* The coordinates (x) were unshifted in update */
1722 /* We will not sum ekinh_old,
1723 * so signal that we still have to do it.
1725 bSumEkinhOld = TRUE;
1730 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1732 /* use the directly determined last velocity, not actually the averaged half steps */
1733 if (bTrotter && ir->eI == IntegrationAlgorithm::VV)
1735 enerd->term[F_EKIN] = last_ekin;
1737 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1739 if (integratorHasConservedEnergyQuantity(ir))
1743 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1747 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1750 /* ######### END PREPARING EDR OUTPUT ########### */
1756 if (fplog && do_log && bDoExpanded)
1758 /* only needed if doing expanded ensemble */
1759 PrintFreeEnergyInfoToFile(fplog,
1761 ir->expandedvals.get(),
1762 ir->bSimTemp ? ir->simtempvals.get() : nullptr,
1763 state_global->dfhist,
1770 energyOutput.addDataAtEnergyStep(bDoDHDL,
1776 ir->expandedvals.get(),
1778 PTCouplingArrays{ state->boxv,
1779 state->nosehoover_xi,
1780 state->nosehoover_vxi,
1782 state->nhpres_vxi },
1794 energyOutput.recordNonEnergyStep();
1797 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1798 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1800 if (doSimulatedAnnealing)
1802 gmx::EnergyOutput::printAnnealingTemperatures(
1803 do_log ? fplog : nullptr, groups, &(ir->opts));
1805 if (do_log || do_ene || do_dr || do_or)
1807 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1811 do_log ? fplog : nullptr,
1817 if (do_log && ir->bDoAwh && awh->hasFepLambdaDimension())
1819 const bool isInitialOutput = false;
1820 printLambdaStateToLog(fplog, state->lambda, isInitialOutput);
1825 pull_print_output(pull_work, step, t);
1828 if (do_per_step(step, ir->nstlog))
1830 if (fflush(fplog) != 0)
1832 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1838 /* Have to do this part _after_ outputting the logfile and the edr file */
1839 /* Gets written into the state at the beginning of next loop*/
1840 state->fep_state = lamnew;
1842 else if (ir->bDoAwh && awh->needForeignEnergyDifferences(step))
1844 state->fep_state = awh->fepLambdaState();
1846 /* Print the remaining wall clock time for the run */
1847 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1851 fprintf(stderr, "\n");
1853 print_time(stderr, walltime_accounting, step, ir, cr);
1856 /* Ion/water position swapping.
1857 * Not done in last step since trajectory writing happens before this call
1858 * in the MD loop and exchanges would be lost anyway. */
1859 bNeedRepartition = FALSE;
1860 if ((ir->eSwapCoords != SwapType::No) && (step > 0) && !bLastStep
1861 && do_per_step(step, ir->swap->nstswap))
1863 bNeedRepartition = do_swapcoords(cr,
1869 as_rvec_array(state->x.data()),
1871 MASTER(cr) && mdrunOptions.verbose,
1874 if (bNeedRepartition && DOMAINDECOMP(cr))
1876 dd_collect_state(cr->dd, state, state_global);
1880 /* Replica exchange */
1884 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1887 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1889 dd_partition_system(fplog,
1910 shouldCheckNumberOfBondedInteractions = true;
1911 upd.setNumAtoms(state->natoms);
1917 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1918 /* With all integrators, except VV, we need to retain the pressure
1919 * at the current step for coupling at the next step.
1921 if ((state->flags & enumValueToBitMask(StateEntry::PressurePrevious))
1922 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1924 /* Store the pressure in t_state for pressure coupling
1925 * at the next MD step.
1927 copy_mat(pres, state->pres_prev);
1930 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1932 if ((membed != nullptr) && (!bLastStep))
1934 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1937 cycles = wallcycle_stop(wcycle, ewcSTEP);
1938 if (DOMAINDECOMP(cr) && wcycle)
1940 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1943 /* increase the MD step number */
1950 fcReportProgress(ir->nsteps + ir->init_step, step);
1954 resetHandler->resetCounters(
1955 step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb, fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1957 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1958 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1960 /* End of main MD loop */
1962 /* Closing TNG files can include compressing data. Therefore it is good to do that
1963 * before stopping the time measurements. */
1964 mdoutf_tng_close(outf);
1966 /* Stop measuring walltime */
1967 walltime_accounting_end_time(walltime_accounting);
1969 if (!thisRankHasDuty(cr, DUTY_PME))
1971 /* Tell the PME only node to finish */
1972 gmx_pme_send_finish(cr);
1977 if (ir->nstcalcenergy > 0)
1979 energyOutput.printEnergyConservation(fplog, ir->simulation_part, EI_MD(ir->eI));
1981 gmx::EnergyOutput::printAnnealingTemperatures(fplog, groups, &(ir->opts));
1982 energyOutput.printAverages(fplog, groups);
1989 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1992 done_shellfc(fplog, shellfc, step_rel);
1994 if (useReplicaExchange && MASTER(cr))
1996 print_replica_exchange_statistics(fplog, repl_ex);
1999 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
2001 global_stat_destroy(gstat);