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 SimulationSignals signals;
204 // Most global communnication stages don't propagate mdrun
205 // signals, and will use this object to achieve that.
206 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
208 if (!mdrunOptions.writeConfout)
210 // This is on by default, and the main known use case for
211 // turning it off is for convenience in benchmarking, which is
212 // something that should not show up in the general user
217 "The -noconfout functionality is deprecated, and may be removed in a "
221 /* md-vv uses averaged full step velocities for T-control
222 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
223 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
224 bTrotter = (EI_VV(ir->eI)
225 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
227 const bool bRerunMD = false;
229 int nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
230 bGStatEveryStep = (nstglobalcomm == 1);
232 const SimulationGroups* groups = &top_global.groups;
234 std::unique_ptr<EssentialDynamics> ed = nullptr;
235 if (opt2bSet("-ei", nfile, fnm))
237 /* Initialize essential dynamics sampling */
238 ed = init_edsam(mdlog,
239 opt2fn_null("-ei", nfile, fnm),
240 opt2fn("-eo", nfile, fnm),
250 else if (observablesHistory->edsamHistory)
253 "The checkpoint is from a run with essential dynamics sampling, "
254 "but the current run did not specify the -ei option. "
255 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
258 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
259 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
260 initialize_lambdas(fplog,
264 ir->simtempvals->temperatures,
265 gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
269 Update upd(*ir, deform);
270 bool doSimulatedAnnealing = false;
272 // TODO: Avoid changing inputrec (#3854)
273 // Simulated annealing updates the reference temperature.
274 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
275 doSimulatedAnnealing = initSimulatedAnnealing(nonConstInputrec, &upd);
277 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
279 const t_fcdata& fcdata = *fr->fcdata;
281 bool simulationsShareState = false;
282 int nstSignalComm = nstglobalcomm;
284 // TODO This implementation of ensemble orientation restraints is nasty because
285 // a user can't just do multi-sim with single-sim orientation restraints.
286 bool usingEnsembleRestraints =
287 (fcdata.disres->nsystems > 1) || ((ms != nullptr) && (fcdata.orires->nr != 0));
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(fr->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->nr == 0,
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 if (ir->efep != FreeEnergyPerturbationType::No)
600 /* Set free energy calculation frequency as the greatest common
601 * denominator of nstdhdl and repl_ex_nst. */
602 nstfep = ir->fepvals->nstdhdl;
605 nstfep = std::gcd(ir->expandedvals->nstexpanded, nstfep);
607 if (useReplicaExchange)
609 nstfep = std::gcd(replExParams.exchangeInterval, nstfep);
613 nstfep = std::gcd(ir->awhParams->nstSampleCoord(), nstfep);
617 /* Be REALLY careful about what flags you set here. You CANNOT assume
618 * this is the first step, since we might be restarting from a checkpoint,
619 * and in that case we should not do any modifications to the state.
621 bStopCM = (ir->comm_mode != ComRemovalAlgorithm::No && !ir->bContinuation);
623 // When restarting from a checkpoint, it can be appropriate to
624 // initialize ekind from quantities in the checkpoint. Otherwise,
625 // compute_globals must initialize ekind before the simulation
626 // starts/restarts. However, only the master rank knows what was
627 // found in the checkpoint file, so we have to communicate in
628 // order to coordinate the restart.
630 // TODO Consider removing this communication if/when checkpoint
631 // reading directly follows .tpr reading, because all ranks can
632 // agree on hasReadEkinState at that time.
633 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
636 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr->mpi_comm_mygroup);
638 if (hasReadEkinState)
640 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
643 unsigned int cglo_flags =
644 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
645 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
647 bSumEkinhOld = FALSE;
649 t_vcm vcm(top_global.groups, *ir);
650 reportComRemovalInfo(fplog, vcm);
652 /* To minimize communication, compute_globals computes the COM velocity
653 * and the kinetic energy for the velocities without COM motion removed.
654 * Thus to get the kinetic energy without the COM contribution, we need
655 * to call compute_globals twice.
657 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
659 unsigned int cglo_flags_iteration = cglo_flags;
660 if (bStopCM && cgloIteration == 0)
662 cglo_flags_iteration |= CGLO_STOPCM;
663 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
665 if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd) && cgloIteration == 0)
667 cglo_flags_iteration |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
669 compute_globals(gstat,
674 makeConstArrayRef(state->x),
675 makeConstArrayRef(state->v),
686 gmx::ArrayRef<real>{},
690 cglo_flags_iteration);
691 if (cglo_flags_iteration & CGLO_STOPCM)
693 /* At initialization, do not pass x with acceleration-correction mode
694 * to avoid (incorrect) correction of the initial coordinates.
696 auto x = (vcm.mode == ComRemovalAlgorithm::LinearAccelerationCorrection)
698 : makeArrayRef(state->x);
699 process_and_stopcm_grp(fplog, &vcm, *md, x, makeArrayRef(state->v));
700 inc_nrnb(nrnb, eNR_STOPCM, md->homenr);
703 if (DOMAINDECOMP(cr))
705 checkNumberOfBondedInteractions(
706 mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
708 if (ir->eI == IntegrationAlgorithm::VVAK)
710 /* a second call to get the half step temperature initialized as well */
711 /* we do the same call as above, but turn the pressure off -- internally to
712 compute_globals, this is recognized as a velocity verlet half-step
713 kinetic energy calculation. This minimized excess variables, but
714 perhaps loses some logic?*/
716 compute_globals(gstat,
721 makeConstArrayRef(state->x),
722 makeConstArrayRef(state->v),
733 gmx::ArrayRef<real>{},
737 cglo_flags & ~CGLO_PRESSURE);
740 /* Calculate the initial half step temperature, and save the ekinh_old */
741 if (startingBehavior == StartingBehavior::NewSimulation)
743 for (i = 0; (i < ir->opts.ngtc); i++)
745 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
749 /* need to make an initiation call to get the Trotter variables set, as well as other constants
750 for non-trotter temperature control */
751 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
755 if (!ir->bContinuation)
757 if (constr && ir->eConstrAlg == ConstraintAlgorithm::Lincs)
760 "RMS relative constraint deviation after constraining: %.2e\n",
763 if (EI_STATE_VELOCITY(ir->eI))
765 real temp = enerd->term[F_TEMP];
766 if (ir->eI != IntegrationAlgorithm::VV)
768 /* Result of Ekin averaged over velocities of -half
769 * and +half step, while we only have -half step here.
773 fprintf(fplog, "Initial temperature: %g K\n", temp);
778 fprintf(stderr, "starting mdrun '%s'\n", *(top_global.name));
781 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
785 sprintf(tbuf, "%s", "infinite");
787 if (ir->init_step > 0)
790 "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
791 gmx_step_str(ir->init_step + ir->nsteps, sbuf),
793 gmx_step_str(ir->init_step, sbuf2),
794 ir->init_step * ir->delta_t);
798 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
800 fprintf(fplog, "\n");
803 walltime_accounting_start_time(walltime_accounting);
804 wallcycle_start(wcycle, WallCycleCounter::Run);
805 print_start(fplog, cr, walltime_accounting, "mdrun");
807 /***********************************************************
811 ************************************************************/
814 /* Skip the first Nose-Hoover integration when we get the state from tpx */
815 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
816 bSumEkinhOld = FALSE;
818 bNeedRepartition = FALSE;
820 step = ir->init_step;
823 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
824 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
825 simulationsShareState,
828 mdrunOptions.reproducible,
830 mdrunOptions.maximumHoursToRun,
835 walltime_accounting);
837 auto checkpointHandler = std::make_unique<CheckpointHandler>(
838 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]),
839 simulationsShareState,
842 mdrunOptions.writeConfout,
843 mdrunOptions.checkpointOptions.period);
845 const bool resetCountersIsLocal = true;
846 auto resetHandler = std::make_unique<ResetHandler>(
847 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
848 !resetCountersIsLocal,
851 mdrunOptions.timingOptions.resetHalfway,
852 mdrunOptions.maximumHoursToRun,
855 walltime_accounting);
857 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
859 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
861 logInitialMultisimStatus(ms, cr, mdlog, simulationsShareState, ir->nsteps, ir->init_step);
864 /* and stop now if we should */
865 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
869 /* Determine if this is a neighbor search step */
870 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
872 if (bPMETune && bNStList)
874 // This has to be here because PME load balancing is called so early.
875 // TODO: Move to after all booleans are defined.
876 if (useGpuForUpdate && !bFirstStep)
878 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
879 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
881 /* PME grid + cut-off optimization with GPUs or PME nodes */
882 pme_loadbal_do(pme_loadbal,
884 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
895 simulationWork.useGpuPmePpCommunication);
898 wallcycle_start(wcycle, WallCycleCounter::Step);
900 bLastStep = (step_rel == ir->nsteps);
901 t = t0 + step * ir->delta_t;
903 // TODO Refactor this, so that nstfep does not need a default value of zero
904 if (ir->efep != FreeEnergyPerturbationType::No || ir->bSimTemp)
906 /* find and set the current lambdas */
907 state->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
909 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
910 bDoFEP = ((ir->efep != FreeEnergyPerturbationType::No) && do_per_step(step, nstfep));
911 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
915 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
916 && do_per_step(step, replExParams.exchangeInterval));
918 if (doSimulatedAnnealing)
920 // TODO: Avoid changing inputrec (#3854)
921 // Simulated annealing updates the reference temperature.
922 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
923 update_annealing_target_temp(nonConstInputrec, t, &upd);
926 /* Stop Center of Mass motion */
927 bStopCM = (ir->comm_mode != ComRemovalAlgorithm::No && do_per_step(step, ir->nstcomm));
929 /* Determine whether or not to do Neighbour Searching */
930 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
932 /* Note that the stopHandler will cause termination at nstglobalcomm
933 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
934 * nstpcouple steps, we have computed the half-step kinetic energy
935 * of the previous step and can always output energies at the last step.
937 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
939 /* do_log triggers energy and virial calculation. Because this leads
940 * to different code paths, forces can be different. Thus for exact
941 * continuation we should avoid extra log output.
942 * Note that the || bLastStep can result in non-exact continuation
943 * beyond the last step. But we don't consider that to be an issue.
945 do_log = (do_per_step(step, ir->nstlog)
946 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
947 do_verbose = mdrunOptions.verbose
948 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
950 if (useGpuForUpdate && !bFirstStep && bNS)
952 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
953 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
954 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
955 // Copy coordinate from the GPU when needed at the search step.
956 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
957 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
958 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
959 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
962 // We only need to calculate virtual velocities if we are writing them in the current step
963 const bool needVirtualVelocitiesThisStep =
965 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep());
967 if (vsite != nullptr)
969 // Virtual sites need to be updated before domain decomposition and forces are calculated
970 wallcycle_start(wcycle, WallCycleCounter::VsiteConstr);
971 // md-vv calculates virtual velocities once it has full-step real velocities
972 vsite->construct(state->x,
975 (!EI_VV(inputrec->eI) && needVirtualVelocitiesThisStep)
976 ? VSiteOperation::PositionsAndVelocities
977 : VSiteOperation::Positions);
978 wallcycle_stop(wcycle, WallCycleCounter::VsiteConstr);
981 if (bNS && !(bFirstStep && ir->bContinuation))
983 bMasterState = FALSE;
984 /* Correct the new box if it is too skewed */
985 if (inputrecDynamicBox(ir))
987 if (correct_box(fplog, step, state->box))
990 // If update is offloaded, it should be informed about the box size change
993 integrator->setPbc(PbcType::Xyz, state->box);
997 if (DOMAINDECOMP(cr) && bMasterState)
999 dd_collect_state(cr->dd, state, state_global);
1002 if (DOMAINDECOMP(cr))
1004 /* Repartition the domain decomposition */
1005 dd_partition_system(fplog,
1025 do_verbose && !bPMETunePrinting);
1026 upd.updateAfterPartition(state->natoms,
1027 md->cFREEZE ? gmx::arrayRefFromArray(md->cFREEZE, md->nr)
1028 : gmx::ArrayRef<const unsigned short>(),
1029 md->cTC ? gmx::arrayRefFromArray(md->cTC, md->nr)
1030 : gmx::ArrayRef<const unsigned short>());
1034 // Allocate or re-size GPU halo exchange object, if necessary
1035 if (bNS && havePPDomainDecomposition(cr) && simulationWork.useGpuHaloExchange)
1037 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
1038 "GPU device manager has to be initialized to use GPU "
1039 "version of halo exchange.");
1040 constructGpuHaloExchange(mdlog, *cr, *fr->deviceStreamManager, wcycle);
1043 if (MASTER(cr) && do_log)
1045 gmx::EnergyOutput::printHeader(
1046 fplog, step, t); /* can we improve the information printed here? */
1049 if (ir->efep != FreeEnergyPerturbationType::No)
1051 update_mdatoms(mdAtoms->mdatoms(), state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
1056 /* We need the kinetic energy at minus the half step for determining
1057 * the full step kinetic energy and possibly for T-coupling.*/
1058 /* This may not be quite working correctly yet . . . . */
1059 int cglo_flags = CGLO_GSTAT | CGLO_TEMPERATURE;
1060 if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
1062 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
1064 compute_globals(gstat,
1069 makeConstArrayRef(state->x),
1070 makeConstArrayRef(state->v),
1081 gmx::ArrayRef<real>{},
1086 if (DOMAINDECOMP(cr))
1088 checkNumberOfBondedInteractions(
1089 mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
1092 clear_mat(force_vir);
1094 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
1096 /* Determine the energy and pressure:
1097 * at nstcalcenergy steps and at energy output steps (set below).
1099 if (EI_VV(ir->eI) && (!bInitStep))
1101 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1102 bCalcVir = bCalcEnerStep
1103 || (ir->epc != PressureCoupling::No
1104 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
1108 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1109 bCalcVir = bCalcEnerStep
1110 || (ir->epc != PressureCoupling::No && do_per_step(step, ir->nstpcouple));
1112 bCalcEner = bCalcEnerStep;
1114 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1116 if (do_ene || do_log || bDoReplEx)
1122 /* Do we need global communication ? */
1123 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
1124 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
1126 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
1127 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
1128 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
1129 if (fr->useMts && !do_per_step(step, ir->nstfout))
1131 force_flags |= GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE;
1136 /* Now is the time to relax the shells */
1137 relax_shell_flexcon(fplog,
1140 mdrunOptions.verbose,
1152 state->x.arrayRefWithPadding(),
1153 state->v.arrayRefWithPadding(),
1168 ddBalanceRegionHandler);
1172 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
1173 is updated (or the AWH update will be performed twice for one step when continuing).
1174 It would be best to call this update function from do_md_trajectory_writing but that
1175 would occur after do_force. One would have to divide the update_awh function into one
1176 function applying the AWH force and one doing the AWH bias update. The update AWH
1177 bias function could then be called after do_md_trajectory_writing (then containing
1178 update_awh_history). The checkpointing will in the future probably moved to the start
1179 of the md loop which will rid of this issue. */
1180 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
1182 awh->updateHistory(state_global->awhHistory.get());
1185 /* The coordinates (x) are shifted (to get whole molecules)
1187 * This is parallellized as well, and does communication too.
1188 * Check comments in sim_util.c
1203 state->x.arrayRefWithPadding(),
1215 ed ? ed->getLegacyED() : nullptr,
1216 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1217 ddBalanceRegionHandler);
1220 // VV integrators do not need the following velocity half step
1221 // if it is the first step after starting from a checkpoint.
1222 // That is, the half step is needed on all other steps, and
1223 // also the first step when starting from a .tpr file.
1226 integrateVVFirstStep(step,
1259 &saved_conserved_quantity,
1269 if (vsite != nullptr && needVirtualVelocitiesThisStep)
1271 // Positions were calculated earlier
1272 wallcycle_start(wcycle, WallCycleCounter::VsiteConstr);
1273 vsite->construct(state->x, state->v, state->box, VSiteOperation::Velocities);
1274 wallcycle_stop(wcycle, WallCycleCounter::VsiteConstr);
1278 /* ######## END FIRST UPDATE STEP ############## */
1279 /* ######## If doing VV, we now have v(dt) ###### */
1282 /* perform extended ensemble sampling in lambda - we don't
1283 actually move to the new state before outputting
1284 statistics, but if performing simulated tempering, we
1285 do update the velocities and the tau_t. */
1286 // TODO: Avoid changing inputrec (#3854)
1287 // Simulated tempering updates the reference temperature.
1288 // Expanded ensemble without simulated tempering does not change the inputrec.
1289 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
1290 lamnew = ExpandedEnsembleDynamics(fplog,
1298 state->v.rvec_array(),
1300 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1303 copy_df_history(state_global->dfhist, state->dfhist);
1307 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1308 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1309 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1310 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1311 || checkpointHandler->isCheckpointingStep()))
1313 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1314 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1316 // Copy velocities if needed for the output/checkpointing.
1317 // NOTE: Copy on the search steps is done at the beginning of the step.
1318 if (useGpuForUpdate && !bNS
1319 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1321 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1322 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1324 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1325 // and update is offloaded hence forces are kept on the GPU for update and have not been
1326 // already transferred in do_force().
1327 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1328 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1329 // prior to GPU update.
1330 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1331 // copy call in do_force(...).
1332 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1333 // on host after the D2H copy in do_force(...).
1334 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1335 && do_per_step(step, ir->nstfout))
1337 stateGpu->copyForcesFromGpu(f.view().force(), AtomLocality::Local);
1338 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1340 /* Now we have the energies and forces corresponding to the
1341 * coordinates at time t. We must output all of this before
1344 do_md_trajectory_writing(fplog,
1361 checkpointHandler->isCheckpointingStep(),
1364 mdrunOptions.writeConfout,
1366 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1367 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x, t);
1369 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1370 if (startingBehavior != StartingBehavior::NewSimulation && bFirstStep
1371 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1373 copy_mat(state->svir_prev, shake_vir);
1374 copy_mat(state->fvir_prev, force_vir);
1377 stopHandler->setSignal();
1378 resetHandler->setSignal(walltime_accounting);
1380 if (bGStat || !PAR(cr))
1382 /* In parallel we only have to check for checkpointing in steps
1383 * where we do global communication,
1384 * otherwise the other nodes don't know.
1386 checkpointHandler->setSignal(walltime_accounting);
1389 /* ######### START SECOND UPDATE STEP ################# */
1391 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1392 controlled in preprocessing */
1394 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1396 gmx_bool bIfRandomize;
1397 bIfRandomize = update_randomize_velocities(ir, step, cr, md, state->v, &upd, constr);
1398 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1399 if (constr && bIfRandomize)
1401 constrain_velocities(constr, do_log, do_ene, step, state, nullptr, false, nullptr);
1404 /* Box is changed in update() when we do pressure coupling,
1405 * but we should still use the old box for energy corrections and when
1406 * writing it to the energy file, so it matches the trajectory files for
1407 * the same timestep above. Make a copy in a separate array.
1409 copy_mat(state->box, lastbox);
1413 if (!useGpuForUpdate)
1415 wallcycle_start(wcycle, WallCycleCounter::Update);
1417 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1420 trotter_update(ir, step, ekind, enerd, state, total_vir, md, &MassQ, trotter_seq, ettTSEQ3);
1421 /* We can only do Berendsen coupling after we have summed
1422 * the kinetic energy or virial. Since the happens
1423 * in global_state after update, we should only do it at
1424 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1429 update_tcouple(step, ir, state, ekind, &MassQ, md);
1430 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1433 /* With leap-frog type integrators we compute the kinetic energy
1434 * at a whole time step as the average of the half-time step kinetic
1435 * energies of two subsequent steps. Therefore we need to compute the
1436 * half step kinetic energy also if we need energies at the next step.
1438 const bool needHalfStepKineticEnergy =
1439 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1441 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1442 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1443 const bool doParrinelloRahman = (ir->epc == PressureCoupling::ParrinelloRahman
1444 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1448 GMX_ASSERT(!useGpuForUpdate, "GPU update is not supported with VVAK integrator.");
1450 integrateVVSecondStep(step,
1486 if (useGpuForUpdate)
1488 if (bNS && (bFirstStep || DOMAINDECOMP(cr)))
1490 integrator->set(stateGpu->getCoordinates(),
1491 stateGpu->getVelocities(),
1492 stateGpu->getForces(),
1496 // Copy data to the GPU after buffers might have being reinitialized
1497 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1498 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1501 if (simulationWork.useGpuPme && !runScheduleWork->simulationWork.useGpuPmePpCommunication
1502 && !thisRankHasDuty(cr, DUTY_PME))
1504 // The PME forces were recieved to the host, so have to be copied
1505 stateGpu->copyForcesToGpu(f.view().force(), AtomLocality::All);
1507 else if (!runScheduleWork->stepWork.useGpuFBufferOps)
1509 // The buffer ops were not offloaded this step, so the forces are on the
1510 // host and have to be copied
1511 stateGpu->copyForcesToGpu(f.view().force(), AtomLocality::Local);
1514 const bool doTemperatureScaling =
1515 (ir->etc != TemperatureCoupling::No
1516 && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1518 // This applies Leap-Frog, LINCS and SETTLE in succession
1519 integrator->integrate(
1520 stateGpu->getForcesReadyOnDeviceEvent(
1521 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1526 doTemperatureScaling,
1529 ir->nstpcouple * ir->delta_t,
1532 // Copy velocities D2H after update if:
1533 // - Globals are computed this step (includes the energy output steps).
1534 // - Temperature is needed for the next step.
1535 if (bGStat || needHalfStepKineticEnergy)
1537 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1538 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1543 /* With multiple time stepping we need to do an additional normal
1544 * update step to obtain the virial, as the actual MTS integration
1545 * using an acceleration where the slow forces are multiplied by mtsFactor.
1546 * Using that acceleration would result in a virial with the slow
1547 * force contribution would be a factor mtsFactor too large.
1549 if (fr->useMts && bCalcVir && constr != nullptr)
1551 upd.update_for_constraint_virial(*ir,
1553 md->havePartiallyFrozenAtoms,
1554 gmx::arrayRefFromArray(md->invmass, md->nr),
1555 gmx::arrayRefFromArray(md->invMassPerDim, md->nr),
1557 f.view().forceWithPadding(),
1560 constrain_coordinates(constr,
1565 upd.xp()->arrayRefWithPadding(),
1571 ArrayRefWithPadding<const RVec> forceCombined =
1572 (fr->useMts && step % ir->mtsLevels[1].stepFactor == 0)
1573 ? f.view().forceMtsCombinedWithPadding()
1574 : f.view().forceWithPadding();
1575 upd.update_coords(*ir,
1578 md->havePartiallyFrozenAtoms,
1579 gmx::arrayRefFromArray(md->ptype, md->nr),
1580 gmx::arrayRefFromArray(md->invmass, md->nr),
1581 gmx::arrayRefFromArray(md->invMassPerDim, md->nr),
1591 wallcycle_stop(wcycle, WallCycleCounter::Update);
1593 constrain_coordinates(constr,
1598 upd.xp()->arrayRefWithPadding(),
1600 bCalcVir && !fr->useMts,
1603 upd.update_sd_second_half(*ir,
1607 gmx::arrayRefFromArray(md->ptype, md->nr),
1608 gmx::arrayRefFromArray(md->invmass, md->nr),
1617 *ir, md->havePartiallyFrozenAtoms, md->homenr, state, wcycle, constr != nullptr);
1620 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1622 updatePrevStepPullCom(pull_work, state);
1625 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1628 /* ############## IF NOT VV, Calculate globals HERE ############ */
1629 /* With Leap-Frog we can skip compute_globals at
1630 * non-communication steps, but we need to calculate
1631 * the kinetic energy one step before communication.
1634 // Organize to do inter-simulation signalling on steps if
1635 // and when algorithms require it.
1636 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1638 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1640 // Copy coordinates when needed to stop the CM motion.
1641 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1643 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1644 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1646 // Since we're already communicating at this step, we
1647 // can propagate intra-simulation signals. Note that
1648 // check_nstglobalcomm has the responsibility for
1649 // choosing the value of nstglobalcomm that is one way
1650 // bGStat becomes true, so we can't get into a
1651 // situation where e.g. checkpointing can't be
1653 bool doIntraSimSignal = true;
1654 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1662 makeConstArrayRef(state->x),
1663 makeConstArrayRef(state->v),
1674 (!EI_VV(ir->eI) && bCalcEner && constr != nullptr) ? constr->rmsdData()
1675 : gmx::ArrayRef<real>{},
1679 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1680 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1681 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1682 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1683 | (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd)
1684 ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1686 if (DOMAINDECOMP(cr))
1688 checkNumberOfBondedInteractions(
1689 mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
1691 if (!EI_VV(ir->eI) && bStopCM)
1693 process_and_stopcm_grp(
1694 fplog, &vcm, *md, makeArrayRef(state->x), makeArrayRef(state->v));
1695 inc_nrnb(nrnb, eNR_STOPCM, md->homenr);
1697 // TODO: The special case of removing CM motion should be dealt more gracefully
1698 if (useGpuForUpdate)
1700 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1701 // Here we block until the H2D copy completes because event sync with the
1702 // force kernels that use the coordinates on the next steps is not implemented
1703 // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1704 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1705 // If the COM removal changed the velocities on the CPU, this has to be accounted for.
1706 if (vcm.mode != ComRemovalAlgorithm::No)
1708 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1715 /* ############# END CALC EKIN AND PRESSURE ################# */
1717 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1718 the virial that should probably be addressed eventually. state->veta has better properies,
1719 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1720 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1722 if (ir->efep != FreeEnergyPerturbationType::No && !EI_VV(ir->eI))
1724 /* Sum up the foreign energy and dK/dl terms for md and sd.
1725 Currently done every step so that dH/dl is correct in the .edr */
1726 accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
1729 update_pcouple_after_coordinates(
1730 fplog, step, ir, md, pres, force_vir, shake_vir, pressureCouplingMu, state, nrnb, upd.deform(), !useGpuForUpdate);
1732 const bool doBerendsenPressureCoupling = (inputrec->epc == PressureCoupling::Berendsen
1733 && do_per_step(step, inputrec->nstpcouple));
1734 const bool doCRescalePressureCoupling = (inputrec->epc == PressureCoupling::CRescale
1735 && do_per_step(step, inputrec->nstpcouple));
1737 && (doBerendsenPressureCoupling || doCRescalePressureCoupling || doParrinelloRahman))
1739 integrator->scaleCoordinates(pressureCouplingMu);
1740 if (doCRescalePressureCoupling)
1742 matrix pressureCouplingInvMu;
1743 gmx::invertBoxMatrix(pressureCouplingMu, pressureCouplingInvMu);
1744 integrator->scaleVelocities(pressureCouplingInvMu);
1746 integrator->setPbc(PbcType::Xyz, state->box);
1749 /* ################# END UPDATE STEP 2 ################# */
1750 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1752 /* The coordinates (x) were unshifted in update */
1755 /* We will not sum ekinh_old,
1756 * so signal that we still have to do it.
1758 bSumEkinhOld = TRUE;
1763 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1765 /* use the directly determined last velocity, not actually the averaged half steps */
1766 if (bTrotter && ir->eI == IntegrationAlgorithm::VV)
1768 enerd->term[F_EKIN] = last_ekin;
1770 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1772 if (integratorHasConservedEnergyQuantity(ir))
1776 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1780 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1783 /* ######### END PREPARING EDR OUTPUT ########### */
1789 if (fplog && do_log && bDoExpanded)
1791 /* only needed if doing expanded ensemble */
1792 PrintFreeEnergyInfoToFile(fplog,
1794 ir->expandedvals.get(),
1795 ir->bSimTemp ? ir->simtempvals.get() : nullptr,
1796 state_global->dfhist,
1803 energyOutput.addDataAtEnergyStep(bDoDHDL,
1809 ir->expandedvals.get(),
1811 PTCouplingArrays{ state->boxv,
1812 state->nosehoover_xi,
1813 state->nosehoover_vxi,
1815 state->nhpres_vxi },
1827 energyOutput.recordNonEnergyStep();
1830 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1831 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1833 if (doSimulatedAnnealing)
1835 gmx::EnergyOutput::printAnnealingTemperatures(
1836 do_log ? fplog : nullptr, groups, &(ir->opts));
1838 if (do_log || do_ene || do_dr || do_or)
1840 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1844 do_log ? fplog : nullptr,
1850 if (do_log && ir->bDoAwh && awh->hasFepLambdaDimension())
1852 const bool isInitialOutput = false;
1853 printLambdaStateToLog(fplog, state->lambda, isInitialOutput);
1858 pull_print_output(pull_work, step, t);
1861 if (do_per_step(step, ir->nstlog))
1863 if (fflush(fplog) != 0)
1865 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1871 /* Have to do this part _after_ outputting the logfile and the edr file */
1872 /* Gets written into the state at the beginning of next loop*/
1873 state->fep_state = lamnew;
1875 else if (ir->bDoAwh && awh->needForeignEnergyDifferences(step))
1877 state->fep_state = awh->fepLambdaState();
1879 /* Print the remaining wall clock time for the run */
1880 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1884 fprintf(stderr, "\n");
1886 print_time(stderr, walltime_accounting, step, ir, cr);
1889 /* Ion/water position swapping.
1890 * Not done in last step since trajectory writing happens before this call
1891 * in the MD loop and exchanges would be lost anyway. */
1892 bNeedRepartition = FALSE;
1893 if ((ir->eSwapCoords != SwapType::No) && (step > 0) && !bLastStep
1894 && do_per_step(step, ir->swap->nstswap))
1896 bNeedRepartition = do_swapcoords(cr,
1902 as_rvec_array(state->x.data()),
1904 MASTER(cr) && mdrunOptions.verbose,
1907 if (bNeedRepartition && DOMAINDECOMP(cr))
1909 dd_collect_state(cr->dd, state, state_global);
1913 /* Replica exchange */
1917 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1920 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1922 dd_partition_system(fplog,
1943 upd.updateAfterPartition(state->natoms,
1944 md->cFREEZE ? gmx::arrayRefFromArray(md->cFREEZE, md->nr)
1945 : gmx::ArrayRef<const unsigned short>(),
1946 md->cTC ? gmx::arrayRefFromArray(md->cTC, md->nr)
1947 : gmx::ArrayRef<const unsigned short>());
1953 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1954 /* With all integrators, except VV, we need to retain the pressure
1955 * at the current step for coupling at the next step.
1957 if ((state->flags & enumValueToBitMask(StateEntry::PressurePrevious))
1958 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1960 /* Store the pressure in t_state for pressure coupling
1961 * at the next MD step.
1963 copy_mat(pres, state->pres_prev);
1966 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1968 if ((membed != nullptr) && (!bLastStep))
1970 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1973 cycles = wallcycle_stop(wcycle, WallCycleCounter::Step);
1974 if (DOMAINDECOMP(cr) && wcycle)
1976 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1979 /* increase the MD step number */
1986 fcReportProgress(ir->nsteps + ir->init_step, step);
1990 resetHandler->resetCounters(
1991 step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb, fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1993 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1994 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1996 /* End of main MD loop */
1998 /* Closing TNG files can include compressing data. Therefore it is good to do that
1999 * before stopping the time measurements. */
2000 mdoutf_tng_close(outf);
2002 /* Stop measuring walltime */
2003 walltime_accounting_end_time(walltime_accounting);
2005 if (!thisRankHasDuty(cr, DUTY_PME))
2007 /* Tell the PME only node to finish */
2008 gmx_pme_send_finish(cr);
2013 if (ir->nstcalcenergy > 0)
2015 energyOutput.printEnergyConservation(fplog, ir->simulation_part, EI_MD(ir->eI));
2017 gmx::EnergyOutput::printAnnealingTemperatures(fplog, groups, &(ir->opts));
2018 energyOutput.printAverages(fplog, groups);
2025 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
2028 done_shellfc(fplog, shellfc, step_rel);
2030 if (useReplicaExchange && MASTER(cr))
2032 print_replica_exchange_statistics(fplog, repl_ex);
2035 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
2037 global_stat_destroy(gstat);