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/commandline/filenm.h"
57 #include "gromacs/domdec/collect.h"
58 #include "gromacs/domdec/dlbtiming.h"
59 #include "gromacs/domdec/domdec.h"
60 #include "gromacs/domdec/domdec_network.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/gpuhaloexchange.h"
63 #include "gromacs/domdec/mdsetup.h"
64 #include "gromacs/domdec/partition.h"
65 #include "gromacs/essentialdynamics/edsam.h"
66 #include "gromacs/ewald/pme_load_balancing.h"
67 #include "gromacs/ewald/pme_pp.h"
68 #include "gromacs/fileio/trxio.h"
69 #include "gromacs/gmxlib/network.h"
70 #include "gromacs/gmxlib/nrnb.h"
71 #include "gromacs/gpu_utils/device_stream_manager.h"
72 #include "gromacs/gpu_utils/gpu_utils.h"
73 #include "gromacs/imd/imd.h"
74 #include "gromacs/listed_forces/listed_forces.h"
75 #include "gromacs/math/functions.h"
76 #include "gromacs/math/invertmatrix.h"
77 #include "gromacs/math/vec.h"
78 #include "gromacs/math/vectypes.h"
79 #include "gromacs/mdlib/checkpointhandler.h"
80 #include "gromacs/mdlib/compute_io.h"
81 #include "gromacs/mdlib/constr.h"
82 #include "gromacs/mdlib/coupling.h"
83 #include "gromacs/mdlib/ebin.h"
84 #include "gromacs/mdlib/enerdata_utils.h"
85 #include "gromacs/mdlib/energyoutput.h"
86 #include "gromacs/mdlib/expanded.h"
87 #include "gromacs/mdlib/force.h"
88 #include "gromacs/mdlib/force_flags.h"
89 #include "gromacs/mdlib/forcerec.h"
90 #include "gromacs/mdlib/freeenergyparameters.h"
91 #include "gromacs/mdlib/md_support.h"
92 #include "gromacs/mdlib/mdatoms.h"
93 #include "gromacs/mdlib/mdoutf.h"
94 #include "gromacs/mdlib/membed.h"
95 #include "gromacs/mdlib/resethandler.h"
96 #include "gromacs/mdlib/sighandler.h"
97 #include "gromacs/mdlib/simulationsignal.h"
98 #include "gromacs/mdlib/stat.h"
99 #include "gromacs/mdlib/stophandler.h"
100 #include "gromacs/mdlib/tgroup.h"
101 #include "gromacs/mdlib/trajectory_writing.h"
102 #include "gromacs/mdlib/update.h"
103 #include "gromacs/mdlib/update_constrain_gpu.h"
104 #include "gromacs/mdlib/update_vv.h"
105 #include "gromacs/mdlib/vcm.h"
106 #include "gromacs/mdlib/vsite.h"
107 #include "gromacs/mdrunutility/handlerestart.h"
108 #include "gromacs/mdrunutility/multisim.h"
109 #include "gromacs/mdrunutility/printtime.h"
110 #include "gromacs/mdtypes/awh_history.h"
111 #include "gromacs/mdtypes/awh_params.h"
112 #include "gromacs/mdtypes/commrec.h"
113 #include "gromacs/mdtypes/df_history.h"
114 #include "gromacs/mdtypes/energyhistory.h"
115 #include "gromacs/mdtypes/fcdata.h"
116 #include "gromacs/mdtypes/forcebuffers.h"
117 #include "gromacs/mdtypes/forcerec.h"
118 #include "gromacs/mdtypes/group.h"
119 #include "gromacs/mdtypes/inputrec.h"
120 #include "gromacs/mdtypes/interaction_const.h"
121 #include "gromacs/mdtypes/md_enums.h"
122 #include "gromacs/mdtypes/mdatom.h"
123 #include "gromacs/mdtypes/mdrunoptions.h"
124 #include "gromacs/mdtypes/multipletimestepping.h"
125 #include "gromacs/mdtypes/observableshistory.h"
126 #include "gromacs/mdtypes/pullhistory.h"
127 #include "gromacs/mdtypes/simulation_workload.h"
128 #include "gromacs/mdtypes/state.h"
129 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
130 #include "gromacs/modularsimulator/energydata.h"
131 #include "gromacs/nbnxm/gpu_data_mgmt.h"
132 #include "gromacs/nbnxm/nbnxm.h"
133 #include "gromacs/pbcutil/pbc.h"
134 #include "gromacs/pulling/output.h"
135 #include "gromacs/pulling/pull.h"
136 #include "gromacs/swap/swapcoords.h"
137 #include "gromacs/timing/wallcycle.h"
138 #include "gromacs/timing/walltime_accounting.h"
139 #include "gromacs/topology/atoms.h"
140 #include "gromacs/topology/idef.h"
141 #include "gromacs/topology/mtop_util.h"
142 #include "gromacs/topology/topology.h"
143 #include "gromacs/trajectory/trajectoryframe.h"
144 #include "gromacs/utility/basedefinitions.h"
145 #include "gromacs/utility/cstringutil.h"
146 #include "gromacs/utility/fatalerror.h"
147 #include "gromacs/utility/logger.h"
148 #include "gromacs/utility/real.h"
149 #include "gromacs/utility/smalloc.h"
151 #include "legacysimulator.h"
152 #include "replicaexchange.h"
155 using gmx::SimulationSignaller;
157 void gmx::LegacySimulator::do_md()
159 // TODO Historically, the EM and MD "integrators" used different
160 // names for the t_inputrec *parameter, but these must have the
161 // same name, now that it's a member of a struct. We use this ir
162 // alias to avoid a large ripple of nearly useless changes.
163 // t_inputrec is being replaced by IMdpOptionsProvider, so this
164 // will go away eventually.
165 const t_inputrec* ir = inputrec;
167 int64_t step, step_rel;
168 double t, t0 = ir->init_t;
169 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
170 gmx_bool bNS = FALSE, bNStList, bStopCM, bFirstStep, bInitStep, bLastStep = FALSE;
171 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
172 gmx_bool do_ene, do_log, do_verbose;
173 gmx_bool bMasterState;
174 unsigned int force_flags;
175 tensor force_vir = { { 0 } }, shake_vir = { { 0 } }, total_vir = { { 0 } }, pres = { { 0 } };
178 matrix pressureCouplingMu, M;
179 gmx_repl_ex_t repl_ex = nullptr;
180 gmx_global_stat_t gstat;
181 gmx_shellfc_t* shellfc;
182 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
185 std::vector<RVec> cbuf;
191 real saved_conserved_quantity = 0;
194 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
196 /* PME load balancing data for GPU kernels */
197 gmx_bool bPMETune = FALSE;
198 gmx_bool bPMETunePrinting = FALSE;
200 bool bInteractiveMDstep = false;
202 /* Domain decomposition could incorrectly miss a bonded
203 interaction, but checking for that requires a global
204 communication stage, which does not otherwise happen in DD
205 code. So we do that alongside the first global energy reduction
206 after a new DD is made. These variables handle whether the
207 check happens, and the result it returns. */
208 bool shouldCheckNumberOfBondedInteractions = false;
209 int totalNumberOfBondedInteractions = -1;
211 SimulationSignals signals;
212 // Most global communnication stages don't propagate mdrun
213 // signals, and will use this object to achieve that.
214 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
216 if (!mdrunOptions.writeConfout)
218 // This is on by default, and the main known use case for
219 // turning it off is for convenience in benchmarking, which is
220 // something that should not show up in the general user
225 "The -noconfout functionality is deprecated, and may be removed in a "
229 /* md-vv uses averaged full step velocities for T-control
230 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
231 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
232 bTrotter = (EI_VV(ir->eI)
233 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
235 const bool bRerunMD = false;
237 int nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
238 bGStatEveryStep = (nstglobalcomm == 1);
240 const SimulationGroups* groups = &top_global->groups;
242 std::unique_ptr<EssentialDynamics> ed = nullptr;
243 if (opt2bSet("-ei", nfile, fnm))
245 /* Initialize essential dynamics sampling */
246 ed = init_edsam(mdlog,
247 opt2fn_null("-ei", nfile, fnm),
248 opt2fn("-eo", nfile, fnm),
258 else if (observablesHistory->edsamHistory)
261 "The checkpoint is from a run with essential dynamics sampling, "
262 "but the current run did not specify the -ei option. "
263 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
266 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
267 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
268 initialize_lambdas(fplog, *ir, MASTER(cr), fep_state, lambda);
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 auto mdatoms = mdAtoms->mdatoms();
368 ForceBuffers f(fr->useMts,
369 ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
370 ? PinningPolicy::PinnedIfSupported
371 : PinningPolicy::CannotBePinned);
372 if (DOMAINDECOMP(cr))
374 stateInstance = std::make_unique<t_state>();
375 state = stateInstance.get();
376 dd_init_local_state(cr->dd, state_global, state);
378 /* Distribute the charge groups over the nodes from the master node */
379 dd_partition_system(fplog,
400 shouldCheckNumberOfBondedInteractions = true;
401 upd.setNumAtoms(state->natoms);
405 state_change_natoms(state_global, state_global->natoms);
406 /* Copy the pointer to the global state */
407 state = state_global;
409 /* Generate and initialize new topology */
410 mdAlgorithmsSetupAtomData(cr, *ir, *top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
412 upd.setNumAtoms(state->natoms);
415 std::unique_ptr<UpdateConstrainGpu> integrator;
417 StatePropagatorDataGpu* stateGpu = fr->stateGpu;
419 // TODO: the assertions below should be handled by UpdateConstraintsBuilder.
422 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
423 || constr->numConstraintsTotal() == 0,
424 "Constraints in domain decomposition are only supported with update "
425 "groups if using GPU update.\n");
426 GMX_RELEASE_ASSERT(ir->eConstrAlg != ConstraintAlgorithm::Shake || constr == nullptr
427 || constr->numConstraintsTotal() == 0,
428 "SHAKE is not supported with GPU update.");
429 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
430 "Either PME or short-ranged non-bonded interaction tasks must run on "
431 "the GPU to use GPU update.\n");
432 GMX_RELEASE_ASSERT(ir->eI == IntegrationAlgorithm::MD,
433 "Only the md integrator is supported with the GPU update.\n");
435 ir->etc != TemperatureCoupling::NoseHoover,
436 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
438 ir->epc == PressureCoupling::No || ir->epc == PressureCoupling::ParrinelloRahman
439 || ir->epc == PressureCoupling::Berendsen || ir->epc == PressureCoupling::CRescale,
440 "Only Parrinello-Rahman, Berendsen, and C-rescale pressure coupling are supported "
441 "with the GPU update.\n");
442 GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
443 "Virtual sites are not supported with the GPU update.\n");
444 GMX_RELEASE_ASSERT(ed == nullptr,
445 "Essential dynamics is not supported with the GPU update.\n");
446 GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(*ir->pull),
447 "Constraints pulling is not supported with the GPU update.\n");
448 GMX_RELEASE_ASSERT(fcdata.orires->nr == 0,
449 "Orientation restraints are not supported with the GPU update.\n");
451 ir->efep == FreeEnergyPerturbationType::No
452 || (!haveFepPerturbedMasses(*top_global) && !havePerturbedConstraints(*top_global)),
453 "Free energy perturbation of masses and constraints are not supported with the GPU "
456 if (constr != nullptr && constr->numConstraintsTotal() > 0)
460 .appendText("Updating coordinates and applying constraints on the GPU.");
464 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
466 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
467 "Device stream manager should be initialized in order to use GPU "
468 "update-constraints.");
470 fr->deviceStreamManager->streamIsValid(gmx::DeviceStreamType::UpdateAndConstraints),
471 "Update stream should be initialized in order to use GPU "
472 "update-constraints.");
473 integrator = std::make_unique<UpdateConstrainGpu>(
477 fr->deviceStreamManager->context(),
478 fr->deviceStreamManager->stream(gmx::DeviceStreamType::UpdateAndConstraints),
479 stateGpu->xUpdatedOnDevice(),
482 integrator->setPbc(PbcType::Xyz, state->box);
485 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
487 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
491 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
494 // NOTE: The global state is no longer used at this point.
495 // But state_global is still used as temporary storage space for writing
496 // the global state to file and potentially for replica exchange.
497 // (Global topology should persist.)
499 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
503 /* Check nstexpanded here, because the grompp check was broken */
504 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
507 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
509 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
514 EnergyData::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
517 preparePrevStepPullCom(
518 ir, pull_work, mdatoms->massT, state, state_global, cr, startingBehavior != StartingBehavior::NewSimulation);
520 // TODO: Remove this by converting AWH into a ForceProvider
521 auto awh = prepareAwhModule(fplog,
526 startingBehavior != StartingBehavior::NewSimulation,
528 opt2fn("-awh", nfile, fnm),
531 if (useReplicaExchange && MASTER(cr))
533 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir, replExParams);
535 /* PME tuning is only supported in the Verlet scheme, with PME for
536 * Coulomb. It is not supported with only LJ PME. */
537 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
538 && ir->cutoff_scheme != CutoffScheme::Group);
540 pme_load_balancing_t* pme_loadbal = nullptr;
544 &pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata, fr->nbv->useGpu());
547 if (!ir->bContinuation)
549 if (state->flags & enumValueToBitMask(StateEntry::V))
551 auto v = makeArrayRef(state->v);
552 /* Set the velocities of vsites, shells and frozen atoms to zero */
553 for (i = 0; i < mdatoms->homenr; i++)
555 if (mdatoms->ptype[i] == eptShell)
559 else if (mdatoms->cFREEZE)
561 for (m = 0; m < DIM; m++)
563 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
574 /* Constrain the initial coordinates and velocities */
575 do_constrain_first(fplog,
580 state->x.arrayRefWithPadding(),
581 state->v.arrayRefWithPadding(),
583 state->lambda[FreeEnergyPerturbationCouplingType::Bonded]);
587 if (ir->efep != FreeEnergyPerturbationType::No)
589 /* Set free energy calculation frequency as the greatest common
590 * denominator of nstdhdl and repl_ex_nst. */
591 nstfep = ir->fepvals->nstdhdl;
594 nstfep = std::gcd(ir->expandedvals->nstexpanded, nstfep);
596 if (useReplicaExchange)
598 nstfep = std::gcd(replExParams.exchangeInterval, nstfep);
602 nstfep = std::gcd(ir->awhParams->nstSampleCoord, nstfep);
606 /* Be REALLY careful about what flags you set here. You CANNOT assume
607 * this is the first step, since we might be restarting from a checkpoint,
608 * and in that case we should not do any modifications to the state.
610 bStopCM = (ir->comm_mode != ComRemovalAlgorithm::No && !ir->bContinuation);
612 // When restarting from a checkpoint, it can be appropriate to
613 // initialize ekind from quantities in the checkpoint. Otherwise,
614 // compute_globals must initialize ekind before the simulation
615 // starts/restarts. However, only the master rank knows what was
616 // found in the checkpoint file, so we have to communicate in
617 // order to coordinate the restart.
619 // TODO Consider removing this communication if/when checkpoint
620 // reading directly follows .tpr reading, because all ranks can
621 // agree on hasReadEkinState at that time.
622 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
625 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr->mpi_comm_mygroup);
627 if (hasReadEkinState)
629 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
632 unsigned int cglo_flags =
633 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
634 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
636 bSumEkinhOld = FALSE;
638 t_vcm vcm(top_global->groups, *ir);
639 reportComRemovalInfo(fplog, vcm);
641 /* To minimize communication, compute_globals computes the COM velocity
642 * and the kinetic energy for the velocities without COM motion removed.
643 * Thus to get the kinetic energy without the COM contribution, we need
644 * to call compute_globals twice.
646 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
648 unsigned int cglo_flags_iteration = cglo_flags;
649 if (bStopCM && cgloIteration == 0)
651 cglo_flags_iteration |= CGLO_STOPCM;
652 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
654 compute_globals(gstat,
659 makeConstArrayRef(state->x),
660 makeConstArrayRef(state->v),
671 gmx::ArrayRef<real>{},
674 &totalNumberOfBondedInteractions,
677 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
679 if (cglo_flags_iteration & CGLO_STOPCM)
681 /* At initialization, do not pass x with acceleration-correction mode
682 * to avoid (incorrect) correction of the initial coordinates.
684 auto x = (vcm.mode == ComRemovalAlgorithm::LinearAccelerationCorrection)
686 : makeArrayRef(state->x);
687 process_and_stopcm_grp(fplog, &vcm, *mdatoms, x, makeArrayRef(state->v));
688 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
691 checkNumberOfBondedInteractions(mdlog,
693 totalNumberOfBondedInteractions,
696 makeConstArrayRef(state->x),
698 &shouldCheckNumberOfBondedInteractions);
699 if (ir->eI == IntegrationAlgorithm::VVAK)
701 /* a second call to get the half step temperature initialized as well */
702 /* we do the same call as above, but turn the pressure off -- internally to
703 compute_globals, this is recognized as a velocity verlet half-step
704 kinetic energy calculation. This minimized excess variables, but
705 perhaps loses some logic?*/
707 compute_globals(gstat,
712 makeConstArrayRef(state->x),
713 makeConstArrayRef(state->v),
724 gmx::ArrayRef<real>{},
729 cglo_flags & ~CGLO_PRESSURE);
732 /* Calculate the initial half step temperature, and save the ekinh_old */
733 if (startingBehavior == StartingBehavior::NewSimulation)
735 for (i = 0; (i < ir->opts.ngtc); i++)
737 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
741 /* need to make an initiation call to get the Trotter variables set, as well as other constants
742 for non-trotter temperature control */
743 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
747 if (!ir->bContinuation)
749 if (constr && ir->eConstrAlg == ConstraintAlgorithm::Lincs)
752 "RMS relative constraint deviation after constraining: %.2e\n",
755 if (EI_STATE_VELOCITY(ir->eI))
757 real temp = enerd->term[F_TEMP];
758 if (ir->eI != IntegrationAlgorithm::VV)
760 /* Result of Ekin averaged over velocities of -half
761 * and +half step, while we only have -half step here.
765 fprintf(fplog, "Initial temperature: %g K\n", temp);
770 fprintf(stderr, "starting mdrun '%s'\n", *(top_global->name));
773 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
777 sprintf(tbuf, "%s", "infinite");
779 if (ir->init_step > 0)
782 "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
783 gmx_step_str(ir->init_step + ir->nsteps, sbuf),
785 gmx_step_str(ir->init_step, sbuf2),
786 ir->init_step * ir->delta_t);
790 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
792 fprintf(fplog, "\n");
795 walltime_accounting_start_time(walltime_accounting);
796 wallcycle_start(wcycle, ewcRUN);
797 print_start(fplog, cr, walltime_accounting, "mdrun");
799 /***********************************************************
803 ************************************************************/
806 /* Skip the first Nose-Hoover integration when we get the state from tpx */
807 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
808 bSumEkinhOld = FALSE;
810 bNeedRepartition = FALSE;
812 step = ir->init_step;
815 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
816 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
817 simulationsShareState,
820 mdrunOptions.reproducible,
822 mdrunOptions.maximumHoursToRun,
827 walltime_accounting);
829 auto checkpointHandler = std::make_unique<CheckpointHandler>(
830 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]),
831 simulationsShareState,
834 mdrunOptions.writeConfout,
835 mdrunOptions.checkpointOptions.period);
837 const bool resetCountersIsLocal = true;
838 auto resetHandler = std::make_unique<ResetHandler>(
839 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
840 !resetCountersIsLocal,
843 mdrunOptions.timingOptions.resetHalfway,
844 mdrunOptions.maximumHoursToRun,
847 walltime_accounting);
849 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
851 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
853 logInitialMultisimStatus(ms, cr, mdlog, simulationsShareState, ir->nsteps, ir->init_step);
856 /* and stop now if we should */
857 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
861 /* Determine if this is a neighbor search step */
862 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
864 if (bPMETune && bNStList)
866 // This has to be here because PME load balancing is called so early.
867 // TODO: Move to after all booleans are defined.
868 if (useGpuForUpdate && !bFirstStep)
870 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
871 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
873 /* PME grid + cut-off optimization with GPUs or PME nodes */
874 pme_loadbal_do(pme_loadbal,
876 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
887 simulationWork.useGpuPmePpCommunication);
890 wallcycle_start(wcycle, ewcSTEP);
892 bLastStep = (step_rel == ir->nsteps);
893 t = t0 + step * ir->delta_t;
895 // TODO Refactor this, so that nstfep does not need a default value of zero
896 if (ir->efep != FreeEnergyPerturbationType::No || ir->bSimTemp)
898 /* find and set the current lambdas */
899 state->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
901 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
902 bDoFEP = ((ir->efep != FreeEnergyPerturbationType::No) && do_per_step(step, nstfep));
903 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
907 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
908 && do_per_step(step, replExParams.exchangeInterval));
910 if (doSimulatedAnnealing)
912 // TODO: Avoid changing inputrec (#3854)
913 // Simulated annealing updates the reference temperature.
914 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
915 update_annealing_target_temp(nonConstInputrec, t, &upd);
918 /* Stop Center of Mass motion */
919 bStopCM = (ir->comm_mode != ComRemovalAlgorithm::No && do_per_step(step, ir->nstcomm));
921 /* Determine whether or not to do Neighbour Searching */
922 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
924 /* Note that the stopHandler will cause termination at nstglobalcomm
925 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
926 * nstpcouple steps, we have computed the half-step kinetic energy
927 * of the previous step and can always output energies at the last step.
929 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
931 /* do_log triggers energy and virial calculation. Because this leads
932 * to different code paths, forces can be different. Thus for exact
933 * continuation we should avoid extra log output.
934 * Note that the || bLastStep can result in non-exact continuation
935 * beyond the last step. But we don't consider that to be an issue.
937 do_log = (do_per_step(step, ir->nstlog)
938 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
939 do_verbose = mdrunOptions.verbose
940 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
942 if (useGpuForUpdate && !bFirstStep && bNS)
944 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
945 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
946 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
947 // Copy coordinate from the GPU when needed at the search step.
948 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
949 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
950 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
951 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
954 // We only need to calculate virtual velocities if we are writing them in the current step
955 const bool needVirtualVelocitiesThisStep =
957 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep());
959 if (vsite != nullptr)
961 // Virtual sites need to be updated before domain decomposition and forces are calculated
962 wallcycle_start(wcycle, ewcVSITECONSTR);
963 // md-vv calculates virtual velocities once it has full-step real velocities
964 vsite->construct(state->x,
967 (!EI_VV(inputrec->eI) && needVirtualVelocitiesThisStep)
968 ? VSiteOperation::PositionsAndVelocities
969 : VSiteOperation::Positions);
970 wallcycle_stop(wcycle, ewcVSITECONSTR);
973 if (bNS && !(bFirstStep && ir->bContinuation))
975 bMasterState = FALSE;
976 /* Correct the new box if it is too skewed */
977 if (inputrecDynamicBox(ir))
979 if (correct_box(fplog, step, state->box))
982 // If update is offloaded, it should be informed about the box size change
985 integrator->setPbc(PbcType::Xyz, state->box);
989 if (DOMAINDECOMP(cr) && bMasterState)
991 dd_collect_state(cr->dd, state, state_global);
994 if (DOMAINDECOMP(cr))
996 /* Repartition the domain decomposition */
997 dd_partition_system(fplog,
1017 do_verbose && !bPMETunePrinting);
1018 shouldCheckNumberOfBondedInteractions = true;
1019 upd.setNumAtoms(state->natoms);
1023 // Allocate or re-size GPU halo exchange object, if necessary
1024 if (bNS && havePPDomainDecomposition(cr) && simulationWork.useGpuHaloExchange)
1026 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
1027 "GPU device manager has to be initialized to use GPU "
1028 "version of halo exchange.");
1029 constructGpuHaloExchange(mdlog, *cr, *fr->deviceStreamManager, wcycle);
1032 if (MASTER(cr) && do_log)
1034 gmx::EnergyOutput::printHeader(
1035 fplog, step, t); /* can we improve the information printed here? */
1038 if (ir->efep != FreeEnergyPerturbationType::No)
1040 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
1045 /* We need the kinetic energy at minus the half step for determining
1046 * the full step kinetic energy and possibly for T-coupling.*/
1047 /* This may not be quite working correctly yet . . . . */
1048 compute_globals(gstat,
1053 makeConstArrayRef(state->x),
1054 makeConstArrayRef(state->v),
1065 gmx::ArrayRef<real>{},
1068 &totalNumberOfBondedInteractions,
1070 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1071 checkNumberOfBondedInteractions(mdlog,
1073 totalNumberOfBondedInteractions,
1076 makeConstArrayRef(state->x),
1078 &shouldCheckNumberOfBondedInteractions);
1080 clear_mat(force_vir);
1082 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
1084 /* Determine the energy and pressure:
1085 * at nstcalcenergy steps and at energy output steps (set below).
1087 if (EI_VV(ir->eI) && (!bInitStep))
1089 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1090 bCalcVir = bCalcEnerStep
1091 || (ir->epc != PressureCoupling::No
1092 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
1096 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1097 bCalcVir = bCalcEnerStep
1098 || (ir->epc != PressureCoupling::No && do_per_step(step, ir->nstpcouple));
1100 bCalcEner = bCalcEnerStep;
1102 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1104 if (do_ene || do_log || bDoReplEx)
1110 /* Do we need global communication ? */
1111 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
1112 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
1114 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
1115 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
1116 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
1117 if (fr->useMts && !do_per_step(step, ir->nstfout))
1119 force_flags |= GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE;
1124 /* Now is the time to relax the shells */
1125 relax_shell_flexcon(fplog,
1128 mdrunOptions.verbose,
1140 state->x.arrayRefWithPadding(),
1141 state->v.arrayRefWithPadding(),
1156 ddBalanceRegionHandler);
1160 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
1161 is updated (or the AWH update will be performed twice for one step when continuing).
1162 It would be best to call this update function from do_md_trajectory_writing but that
1163 would occur after do_force. One would have to divide the update_awh function into one
1164 function applying the AWH force and one doing the AWH bias update. The update AWH
1165 bias function could then be called after do_md_trajectory_writing (then containing
1166 update_awh_history). The checkpointing will in the future probably moved to the start
1167 of the md loop which will rid of this issue. */
1168 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
1170 awh->updateHistory(state_global->awhHistory.get());
1173 /* The coordinates (x) are shifted (to get whole molecules)
1175 * This is parallellized as well, and does communication too.
1176 * Check comments in sim_util.c
1191 state->x.arrayRefWithPadding(),
1203 ed ? ed->getLegacyED() : nullptr,
1204 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1205 ddBalanceRegionHandler);
1208 // VV integrators do not need the following velocity half step
1209 // if it is the first step after starting from a checkpoint.
1210 // That is, the half step is needed on all other steps, and
1211 // also the first step when starting from a .tpr file.
1214 integrateVVFirstStep(step,
1247 &shouldCheckNumberOfBondedInteractions,
1248 &saved_conserved_quantity,
1258 if (vsite != nullptr && needVirtualVelocitiesThisStep)
1260 // Positions were calculated earlier
1261 wallcycle_start(wcycle, ewcVSITECONSTR);
1262 vsite->construct(state->x, state->v, state->box, VSiteOperation::Velocities);
1263 wallcycle_stop(wcycle, ewcVSITECONSTR);
1267 /* ######## END FIRST UPDATE STEP ############## */
1268 /* ######## If doing VV, we now have v(dt) ###### */
1271 /* perform extended ensemble sampling in lambda - we don't
1272 actually move to the new state before outputting
1273 statistics, but if performing simulated tempering, we
1274 do update the velocities and the tau_t. */
1275 // TODO: Avoid changing inputrec (#3854)
1276 // Simulated tempering updates the reference temperature.
1277 // Expanded ensemble without simulated tempering does not change the inputrec.
1278 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
1279 lamnew = ExpandedEnsembleDynamics(fplog,
1287 state->v.rvec_array(),
1289 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1292 copy_df_history(state_global->dfhist, state->dfhist);
1296 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1297 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1298 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1299 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1300 || checkpointHandler->isCheckpointingStep()))
1302 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1303 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1305 // Copy velocities if needed for the output/checkpointing.
1306 // NOTE: Copy on the search steps is done at the beginning of the step.
1307 if (useGpuForUpdate && !bNS
1308 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1310 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1311 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1313 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1314 // and update is offloaded hence forces are kept on the GPU for update and have not been
1315 // already transferred in do_force().
1316 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1317 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1318 // prior to GPU update.
1319 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1320 // copy call in do_force(...).
1321 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1322 // on host after the D2H copy in do_force(...).
1323 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1324 && do_per_step(step, ir->nstfout))
1326 stateGpu->copyForcesFromGpu(f.view().force(), AtomLocality::Local);
1327 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1329 /* Now we have the energies and forces corresponding to the
1330 * coordinates at time t. We must output all of this before
1333 do_md_trajectory_writing(fplog,
1350 checkpointHandler->isCheckpointingStep(),
1353 mdrunOptions.writeConfout,
1355 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1356 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1358 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1359 if (startingBehavior != StartingBehavior::NewSimulation && bFirstStep
1360 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1362 copy_mat(state->svir_prev, shake_vir);
1363 copy_mat(state->fvir_prev, force_vir);
1366 stopHandler->setSignal();
1367 resetHandler->setSignal(walltime_accounting);
1369 if (bGStat || !PAR(cr))
1371 /* In parallel we only have to check for checkpointing in steps
1372 * where we do global communication,
1373 * otherwise the other nodes don't know.
1375 checkpointHandler->setSignal(walltime_accounting);
1378 /* ######### START SECOND UPDATE STEP ################# */
1380 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1381 controlled in preprocessing */
1383 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1385 gmx_bool bIfRandomize;
1386 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1387 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1388 if (constr && bIfRandomize)
1390 constrain_velocities(constr, do_log, do_ene, step, state, nullptr, false, nullptr);
1393 /* Box is changed in update() when we do pressure coupling,
1394 * but we should still use the old box for energy corrections and when
1395 * writing it to the energy file, so it matches the trajectory files for
1396 * the same timestep above. Make a copy in a separate array.
1398 copy_mat(state->box, lastbox);
1402 if (!useGpuForUpdate)
1404 wallcycle_start(wcycle, ewcUPDATE);
1406 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1409 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1410 /* We can only do Berendsen coupling after we have summed
1411 * the kinetic energy or virial. Since the happens
1412 * in global_state after update, we should only do it at
1413 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1418 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1419 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1422 /* With leap-frog type integrators we compute the kinetic energy
1423 * at a whole time step as the average of the half-time step kinetic
1424 * energies of two subsequent steps. Therefore we need to compute the
1425 * half step kinetic energy also if we need energies at the next step.
1427 const bool needHalfStepKineticEnergy =
1428 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1430 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1431 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1432 const bool doParrinelloRahman = (ir->epc == PressureCoupling::ParrinelloRahman
1433 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1437 GMX_ASSERT(!useGpuForUpdate, "GPU update is not supported with VVAK integrator.");
1439 integrateVVSecondStep(step,
1475 if (useGpuForUpdate)
1478 wallcycle_stop(wcycle, ewcUPDATE);
1480 if (bNS && (bFirstStep || DOMAINDECOMP(cr)))
1482 integrator->set(stateGpu->getCoordinates(),
1483 stateGpu->getVelocities(),
1484 stateGpu->getForces(),
1488 // Copy data to the GPU after buffers might have being reinitialized
1489 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1490 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1493 if (simulationWork.useGpuPme && !runScheduleWork->simulationWork.useGpuPmePpCommunication
1494 && !thisRankHasDuty(cr, DUTY_PME))
1496 // The PME forces were recieved to the host, so have to be copied
1497 stateGpu->copyForcesToGpu(f.view().force(), AtomLocality::All);
1499 else if (!runScheduleWork->stepWork.useGpuFBufferOps)
1501 // The buffer ops were not offloaded this step, so the forces are on the
1502 // host and have to be copied
1503 stateGpu->copyForcesToGpu(f.view().force(), AtomLocality::Local);
1506 const bool doTemperatureScaling =
1507 (ir->etc != TemperatureCoupling::No
1508 && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1510 // This applies Leap-Frog, LINCS and SETTLE in succession
1511 integrator->integrate(
1512 stateGpu->getForcesReadyOnDeviceEvent(
1513 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1518 doTemperatureScaling,
1521 ir->nstpcouple * ir->delta_t,
1524 // Copy velocities D2H after update if:
1525 // - Globals are computed this step (includes the energy output steps).
1526 // - Temperature is needed for the next step.
1527 if (bGStat || needHalfStepKineticEnergy)
1529 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1530 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1535 /* With multiple time stepping we need to do an additional normal
1536 * update step to obtain the virial, as the actual MTS integration
1537 * using an acceleration where the slow forces are multiplied by mtsFactor.
1538 * Using that acceleration would result in a virial with the slow
1539 * force contribution would be a factor mtsFactor too large.
1541 if (fr->useMts && bCalcVir && constr != nullptr)
1543 upd.update_for_constraint_virial(
1544 *ir, *mdatoms, *state, f.view().forceWithPadding(), *ekind);
1546 constrain_coordinates(constr,
1551 upd.xp()->arrayRefWithPadding(),
1557 ArrayRefWithPadding<const RVec> forceCombined =
1558 (fr->useMts && step % ir->mtsLevels[1].stepFactor == 0)
1559 ? f.view().forceMtsCombinedWithPadding()
1560 : f.view().forceWithPadding();
1562 *ir, step, mdatoms, state, forceCombined, fcdata, ekind, M, etrtPOSITION, cr, constr != nullptr);
1564 wallcycle_stop(wcycle, ewcUPDATE);
1566 constrain_coordinates(constr,
1571 upd.xp()->arrayRefWithPadding(),
1573 bCalcVir && !fr->useMts,
1576 upd.update_sd_second_half(
1577 *ir, step, &dvdl_constr, mdatoms, state, cr, nrnb, wcycle, constr, do_log, do_ene);
1578 upd.finish_update(*ir, mdatoms, state, wcycle, constr != nullptr);
1581 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1583 updatePrevStepPullCom(pull_work, state);
1586 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1589 /* ############## IF NOT VV, Calculate globals HERE ############ */
1590 /* With Leap-Frog we can skip compute_globals at
1591 * non-communication steps, but we need to calculate
1592 * the kinetic energy one step before communication.
1595 // Organize to do inter-simulation signalling on steps if
1596 // and when algorithms require it.
1597 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1599 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1601 // Copy coordinates when needed to stop the CM motion.
1602 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1604 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1605 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1607 // Since we're already communicating at this step, we
1608 // can propagate intra-simulation signals. Note that
1609 // check_nstglobalcomm has the responsibility for
1610 // choosing the value of nstglobalcomm that is one way
1611 // bGStat becomes true, so we can't get into a
1612 // situation where e.g. checkpointing can't be
1614 bool doIntraSimSignal = true;
1615 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1623 makeConstArrayRef(state->x),
1624 makeConstArrayRef(state->v),
1635 (!EI_VV(ir->eI) && bCalcEner && constr != nullptr) ? constr->rmsdData()
1636 : gmx::ArrayRef<real>{},
1639 &totalNumberOfBondedInteractions,
1641 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1642 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1643 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1644 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1645 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1647 checkNumberOfBondedInteractions(mdlog,
1649 totalNumberOfBondedInteractions,
1652 makeConstArrayRef(state->x),
1654 &shouldCheckNumberOfBondedInteractions);
1655 if (!EI_VV(ir->eI) && bStopCM)
1657 process_and_stopcm_grp(
1658 fplog, &vcm, *mdatoms, makeArrayRef(state->x), makeArrayRef(state->v));
1659 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1661 // TODO: The special case of removing CM motion should be dealt more gracefully
1662 if (useGpuForUpdate)
1664 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1665 // Here we block until the H2D copy completes because event sync with the
1666 // force kernels that use the coordinates on the next steps is not implemented
1667 // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1668 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1669 // If the COM removal changed the velocities on the CPU, this has to be accounted for.
1670 if (vcm.mode != ComRemovalAlgorithm::No)
1672 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1679 /* ############# END CALC EKIN AND PRESSURE ################# */
1681 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1682 the virial that should probably be addressed eventually. state->veta has better properies,
1683 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1684 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1686 if (ir->efep != FreeEnergyPerturbationType::No && !EI_VV(ir->eI))
1688 /* Sum up the foreign energy and dK/dl terms for md and sd.
1689 Currently done every step so that dH/dl is correct in the .edr */
1690 accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
1693 update_pcouple_after_coordinates(
1694 fplog, step, ir, mdatoms, pres, force_vir, shake_vir, pressureCouplingMu, state, nrnb, upd.deform(), !useGpuForUpdate);
1696 const bool doBerendsenPressureCoupling = (inputrec->epc == PressureCoupling::Berendsen
1697 && do_per_step(step, inputrec->nstpcouple));
1698 const bool doCRescalePressureCoupling = (inputrec->epc == PressureCoupling::CRescale
1699 && do_per_step(step, inputrec->nstpcouple));
1701 && (doBerendsenPressureCoupling || doCRescalePressureCoupling || doParrinelloRahman))
1703 integrator->scaleCoordinates(pressureCouplingMu);
1704 if (doCRescalePressureCoupling)
1706 matrix pressureCouplingInvMu;
1707 gmx::invertBoxMatrix(pressureCouplingMu, pressureCouplingInvMu);
1708 integrator->scaleVelocities(pressureCouplingInvMu);
1710 integrator->setPbc(PbcType::Xyz, state->box);
1713 /* ################# END UPDATE STEP 2 ################# */
1714 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1716 /* The coordinates (x) were unshifted in update */
1719 /* We will not sum ekinh_old,
1720 * so signal that we still have to do it.
1722 bSumEkinhOld = TRUE;
1727 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1729 /* use the directly determined last velocity, not actually the averaged half steps */
1730 if (bTrotter && ir->eI == IntegrationAlgorithm::VV)
1732 enerd->term[F_EKIN] = last_ekin;
1734 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1736 if (integratorHasConservedEnergyQuantity(ir))
1740 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1744 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1747 /* ######### END PREPARING EDR OUTPUT ########### */
1753 if (fplog && do_log && bDoExpanded)
1755 /* only needed if doing expanded ensemble */
1756 PrintFreeEnergyInfoToFile(fplog,
1758 ir->expandedvals.get(),
1759 ir->bSimTemp ? ir->simtempvals.get() : nullptr,
1760 state_global->dfhist,
1767 energyOutput.addDataAtEnergyStep(bDoDHDL,
1773 ir->expandedvals.get(),
1775 PTCouplingArrays{ state->boxv,
1776 state->nosehoover_xi,
1777 state->nosehoover_vxi,
1779 state->nhpres_vxi },
1791 energyOutput.recordNonEnergyStep();
1794 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1795 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1797 if (doSimulatedAnnealing)
1799 gmx::EnergyOutput::printAnnealingTemperatures(
1800 do_log ? fplog : nullptr, groups, &(ir->opts));
1802 if (do_log || do_ene || do_dr || do_or)
1804 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1808 do_log ? fplog : nullptr,
1814 if (do_log && ir->bDoAwh && awh->hasFepLambdaDimension())
1816 const bool isInitialOutput = false;
1817 printLambdaStateToLog(fplog, state->lambda, isInitialOutput);
1822 pull_print_output(pull_work, step, t);
1825 if (do_per_step(step, ir->nstlog))
1827 if (fflush(fplog) != 0)
1829 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1835 /* Have to do this part _after_ outputting the logfile and the edr file */
1836 /* Gets written into the state at the beginning of next loop*/
1837 state->fep_state = lamnew;
1839 else if (ir->bDoAwh && awh->needForeignEnergyDifferences(step))
1841 state->fep_state = awh->fepLambdaState();
1843 /* Print the remaining wall clock time for the run */
1844 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1848 fprintf(stderr, "\n");
1850 print_time(stderr, walltime_accounting, step, ir, cr);
1853 /* Ion/water position swapping.
1854 * Not done in last step since trajectory writing happens before this call
1855 * in the MD loop and exchanges would be lost anyway. */
1856 bNeedRepartition = FALSE;
1857 if ((ir->eSwapCoords != SwapType::No) && (step > 0) && !bLastStep
1858 && do_per_step(step, ir->swap->nstswap))
1860 bNeedRepartition = do_swapcoords(cr,
1866 as_rvec_array(state->x.data()),
1868 MASTER(cr) && mdrunOptions.verbose,
1871 if (bNeedRepartition && DOMAINDECOMP(cr))
1873 dd_collect_state(cr->dd, state, state_global);
1877 /* Replica exchange */
1881 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1884 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1886 dd_partition_system(fplog,
1907 shouldCheckNumberOfBondedInteractions = true;
1908 upd.setNumAtoms(state->natoms);
1914 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1915 /* With all integrators, except VV, we need to retain the pressure
1916 * at the current step for coupling at the next step.
1918 if ((state->flags & enumValueToBitMask(StateEntry::PressurePrevious))
1919 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1921 /* Store the pressure in t_state for pressure coupling
1922 * at the next MD step.
1924 copy_mat(pres, state->pres_prev);
1927 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1929 if ((membed != nullptr) && (!bLastStep))
1931 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1934 cycles = wallcycle_stop(wcycle, ewcSTEP);
1935 if (DOMAINDECOMP(cr) && wcycle)
1937 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1940 /* increase the MD step number */
1947 fcReportProgress(ir->nsteps + ir->init_step, step);
1951 resetHandler->resetCounters(
1952 step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb, fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1954 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1955 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1957 /* End of main MD loop */
1959 /* Closing TNG files can include compressing data. Therefore it is good to do that
1960 * before stopping the time measurements. */
1961 mdoutf_tng_close(outf);
1963 /* Stop measuring walltime */
1964 walltime_accounting_end_time(walltime_accounting);
1966 if (!thisRankHasDuty(cr, DUTY_PME))
1968 /* Tell the PME only node to finish */
1969 gmx_pme_send_finish(cr);
1974 if (ir->nstcalcenergy > 0)
1976 energyOutput.printEnergyConservation(fplog, ir->simulation_part, EI_MD(ir->eI));
1978 gmx::EnergyOutput::printAnnealingTemperatures(fplog, groups, &(ir->opts));
1979 energyOutput.printAverages(fplog, groups);
1986 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1989 done_shellfc(fplog, shellfc, step_rel);
1991 if (useReplicaExchange && MASTER(cr))
1993 print_replica_exchange_statistics(fplog, repl_ex);
1996 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1998 global_stat_destroy(gstat);