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, 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/vec.h"
77 #include "gromacs/math/vectypes.h"
78 #include "gromacs/mdlib/checkpointhandler.h"
79 #include "gromacs/mdlib/compute_io.h"
80 #include "gromacs/mdlib/constr.h"
81 #include "gromacs/mdlib/coupling.h"
82 #include "gromacs/mdlib/ebin.h"
83 #include "gromacs/mdlib/enerdata_utils.h"
84 #include "gromacs/mdlib/energyoutput.h"
85 #include "gromacs/mdlib/expanded.h"
86 #include "gromacs/mdlib/force.h"
87 #include "gromacs/mdlib/force_flags.h"
88 #include "gromacs/mdlib/forcerec.h"
89 #include "gromacs/mdlib/freeenergyparameters.h"
90 #include "gromacs/mdlib/md_support.h"
91 #include "gromacs/mdlib/mdatoms.h"
92 #include "gromacs/mdlib/mdoutf.h"
93 #include "gromacs/mdlib/membed.h"
94 #include "gromacs/mdlib/resethandler.h"
95 #include "gromacs/mdlib/sighandler.h"
96 #include "gromacs/mdlib/simulationsignal.h"
97 #include "gromacs/mdlib/stat.h"
98 #include "gromacs/mdlib/stophandler.h"
99 #include "gromacs/mdlib/tgroup.h"
100 #include "gromacs/mdlib/trajectory_writing.h"
101 #include "gromacs/mdlib/update.h"
102 #include "gromacs/mdlib/update_constrain_gpu.h"
103 #include "gromacs/mdlib/vcm.h"
104 #include "gromacs/mdlib/vsite.h"
105 #include "gromacs/mdrunutility/handlerestart.h"
106 #include "gromacs/mdrunutility/multisim.h"
107 #include "gromacs/mdrunutility/printtime.h"
108 #include "gromacs/mdtypes/awh_history.h"
109 #include "gromacs/mdtypes/awh_params.h"
110 #include "gromacs/mdtypes/commrec.h"
111 #include "gromacs/mdtypes/df_history.h"
112 #include "gromacs/mdtypes/energyhistory.h"
113 #include "gromacs/mdtypes/fcdata.h"
114 #include "gromacs/mdtypes/forcebuffers.h"
115 #include "gromacs/mdtypes/forcerec.h"
116 #include "gromacs/mdtypes/group.h"
117 #include "gromacs/mdtypes/inputrec.h"
118 #include "gromacs/mdtypes/interaction_const.h"
119 #include "gromacs/mdtypes/md_enums.h"
120 #include "gromacs/mdtypes/mdatom.h"
121 #include "gromacs/mdtypes/mdrunoptions.h"
122 #include "gromacs/mdtypes/observableshistory.h"
123 #include "gromacs/mdtypes/pullhistory.h"
124 #include "gromacs/mdtypes/simulation_workload.h"
125 #include "gromacs/mdtypes/state.h"
126 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
127 #include "gromacs/modularsimulator/energydata.h"
128 #include "gromacs/nbnxm/gpu_data_mgmt.h"
129 #include "gromacs/nbnxm/nbnxm.h"
130 #include "gromacs/pbcutil/pbc.h"
131 #include "gromacs/pulling/output.h"
132 #include "gromacs/pulling/pull.h"
133 #include "gromacs/swap/swapcoords.h"
134 #include "gromacs/timing/wallcycle.h"
135 #include "gromacs/timing/walltime_accounting.h"
136 #include "gromacs/topology/atoms.h"
137 #include "gromacs/topology/idef.h"
138 #include "gromacs/topology/mtop_util.h"
139 #include "gromacs/topology/topology.h"
140 #include "gromacs/trajectory/trajectoryframe.h"
141 #include "gromacs/utility/basedefinitions.h"
142 #include "gromacs/utility/cstringutil.h"
143 #include "gromacs/utility/fatalerror.h"
144 #include "gromacs/utility/logger.h"
145 #include "gromacs/utility/real.h"
146 #include "gromacs/utility/smalloc.h"
148 #include "legacysimulator.h"
149 #include "replicaexchange.h"
153 # include "corewrap.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 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;
183 gmx_bool bTemp, bPres, bTrotter;
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, opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm), top_global,
247 ir, cr, constr, state_global, observablesHistory, oenv, startingBehavior);
249 else if (observablesHistory->edsamHistory)
252 "The checkpoint is from a run with essential dynamics sampling, "
253 "but the current run did not specify the -ei option. "
254 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
257 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
258 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
259 initialize_lambdas(fplog, *ir, MASTER(cr), fep_state, lambda);
260 Update upd(*ir, deform);
261 const bool doSimulatedAnnealing = initSimulatedAnnealing(ir, &upd);
262 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
264 const t_fcdata& fcdata = *fr->fcdata;
266 bool simulationsShareState = false;
267 int nstSignalComm = nstglobalcomm;
269 // TODO This implementation of ensemble orientation restraints is nasty because
270 // a user can't just do multi-sim with single-sim orientation restraints.
271 bool usingEnsembleRestraints =
272 (fcdata.disres->nsystems > 1) || ((ms != nullptr) && (fcdata.orires->nr != 0));
273 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
275 // Replica exchange, ensemble restraints and AWH need all
276 // simulations to remain synchronized, so they need
277 // checkpoints and stop conditions to act on the same step, so
278 // the propagation of such signals must take place between
279 // simulations, not just within simulations.
280 // TODO: Make algorithm initializers set these flags.
281 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
283 if (simulationsShareState)
285 // Inter-simulation signal communication does not need to happen
286 // often, so we use a minimum of 200 steps to reduce overhead.
287 const int c_minimumInterSimulationSignallingInterval = 200;
288 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
293 if (startingBehavior != StartingBehavior::RestartWithAppending)
295 pleaseCiteCouplingAlgorithms(fplog, *ir);
298 init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, ir,
299 top_global, oenv, wcycle, startingBehavior, simulationsShareState, ms);
300 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
301 mdoutf_get_fp_dhdl(outf), false, startingBehavior, mdModulesNotifier);
303 gstat = global_stat_init(ir);
305 /* Check for polarizable models and flexible constraints */
306 shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
307 ir->nstcalcenergy, DOMAINDECOMP(cr));
310 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
311 if ((io > 2000) && MASTER(cr))
313 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
317 // Local state only becomes valid now.
318 std::unique_ptr<t_state> stateInstance;
322 gmx_localtop_t top(top_global->ffparams);
324 auto mdatoms = mdAtoms->mdatoms();
326 const auto& simulationWork = runScheduleWork->simulationWork;
327 const bool useGpuForPme = simulationWork.useGpuPme;
328 const bool useGpuForNonbonded = simulationWork.useGpuNonbonded;
329 const bool useGpuForBufferOps = simulationWork.useGpuBufferOps;
330 const bool useGpuForUpdate = simulationWork.useGpuUpdate;
332 ForceBuffers f(((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
333 ? PinningPolicy::PinnedIfSupported
334 : PinningPolicy::CannotBePinned);
335 if (DOMAINDECOMP(cr))
337 stateInstance = std::make_unique<t_state>();
338 state = stateInstance.get();
339 dd_init_local_state(cr->dd, state_global, state);
341 /* Distribute the charge groups over the nodes from the master node */
342 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
343 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
344 nrnb, nullptr, FALSE);
345 shouldCheckNumberOfBondedInteractions = true;
346 upd.setNumAtoms(state->natoms);
350 state_change_natoms(state_global, state_global->natoms);
351 /* Copy the pointer to the global state */
352 state = state_global;
354 /* Generate and initialize new topology */
355 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
357 upd.setNumAtoms(state->natoms);
360 std::unique_ptr<UpdateConstrainGpu> integrator;
362 StatePropagatorDataGpu* stateGpu = fr->stateGpu;
364 // TODO: the assertions below should be handled by UpdateConstraintsBuilder.
367 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
368 || constr->numConstraintsTotal() == 0,
369 "Constraints in domain decomposition are only supported with update "
370 "groups if using GPU update.\n");
371 GMX_RELEASE_ASSERT(ir->eConstrAlg != econtSHAKE || constr == nullptr
372 || constr->numConstraintsTotal() == 0,
373 "SHAKE is not supported with GPU update.");
374 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
375 "Either PME or short-ranged non-bonded interaction tasks must run on "
376 "the GPU to use GPU update.\n");
377 GMX_RELEASE_ASSERT(ir->eI == eiMD,
378 "Only the md integrator is supported with the GPU update.\n");
380 ir->etc != etcNOSEHOOVER,
381 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
382 GMX_RELEASE_ASSERT(ir->epc == epcNO || ir->epc == epcPARRINELLORAHMAN || ir->epc == epcBERENDSEN,
383 "Only Parrinello-Rahman and Berendsen pressure coupling are supported "
384 "with the GPU update.\n");
385 GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
386 "Virtual sites are not supported with the GPU update.\n");
387 GMX_RELEASE_ASSERT(ed == nullptr,
388 "Essential dynamics is not supported with the GPU update.\n");
389 GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(ir->pull),
390 "Constraints pulling is not supported with the GPU update.\n");
391 GMX_RELEASE_ASSERT(fcdata.orires->nr == 0,
392 "Orientation restraints are not supported with the GPU update.\n");
395 || (!haveFreeEnergyType(*ir, efptBONDED) && !haveFreeEnergyType(*ir, efptMASS)),
396 "Free energy perturbation of masses and constraints are not supported with the GPU "
399 if (constr != nullptr && constr->numConstraintsTotal() > 0)
403 .appendText("Updating coordinates and applying constraints on the GPU.");
407 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
409 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
410 "Device stream manager should be initialized in order to use GPU "
411 "update-constraints.");
413 fr->deviceStreamManager->streamIsValid(gmx::DeviceStreamType::UpdateAndConstraints),
414 "Update stream should be initialized in order to use GPU "
415 "update-constraints.");
416 integrator = std::make_unique<UpdateConstrainGpu>(
417 *ir, *top_global, fr->deviceStreamManager->context(),
418 fr->deviceStreamManager->stream(gmx::DeviceStreamType::UpdateAndConstraints),
419 stateGpu->xUpdatedOnDevice());
421 integrator->setPbc(PbcType::Xyz, state->box);
424 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
426 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
430 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
433 // NOTE: The global state is no longer used at this point.
434 // But state_global is still used as temporary storage space for writing
435 // the global state to file and potentially for replica exchange.
436 // (Global topology should persist.)
438 update_mdatoms(mdatoms, state->lambda[efptMASS]);
442 /* Check nstexpanded here, because the grompp check was broken */
443 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
446 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
448 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
453 EnergyData::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
456 preparePrevStepPullCom(ir, pull_work, mdatoms->massT, state, state_global, cr,
457 startingBehavior != StartingBehavior::NewSimulation);
459 // TODO: Remove this by converting AWH into a ForceProvider
460 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms,
461 startingBehavior != StartingBehavior::NewSimulation,
462 shellfc != nullptr, opt2fn("-awh", nfile, fnm), pull_work);
464 if (useReplicaExchange && MASTER(cr))
466 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir, replExParams);
468 /* PME tuning is only supported in the Verlet scheme, with PME for
469 * Coulomb. It is not supported with only LJ PME. */
470 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
471 && ir->cutoff_scheme != ecutsGROUP);
473 pme_load_balancing_t* pme_loadbal = nullptr;
476 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata,
480 if (!ir->bContinuation)
482 if (state->flags & (1U << estV))
484 auto v = makeArrayRef(state->v);
485 /* Set the velocities of vsites, shells and frozen atoms to zero */
486 for (i = 0; i < mdatoms->homenr; i++)
488 if (mdatoms->ptype[i] == eptVSite || mdatoms->ptype[i] == eptShell)
492 else if (mdatoms->cFREEZE)
494 for (m = 0; m < DIM; m++)
496 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
507 /* Constrain the initial coordinates and velocities */
508 do_constrain_first(fplog, constr, ir, mdatoms->nr, mdatoms->homenr,
509 state->x.arrayRefWithPadding(), state->v.arrayRefWithPadding(),
510 state->box, state->lambda[efptBONDED]);
514 /* Construct the virtual sites for the initial configuration */
515 vsite->construct(state->x, ir->delta_t, {}, state->box);
519 if (ir->efep != efepNO)
521 /* Set free energy calculation frequency as the greatest common
522 * denominator of nstdhdl and repl_ex_nst. */
523 nstfep = ir->fepvals->nstdhdl;
526 nstfep = std::gcd(ir->expandedvals->nstexpanded, nstfep);
528 if (useReplicaExchange)
530 nstfep = std::gcd(replExParams.exchangeInterval, nstfep);
534 nstfep = std::gcd(ir->awhParams->nstSampleCoord, nstfep);
538 /* Be REALLY careful about what flags you set here. You CANNOT assume
539 * this is the first step, since we might be restarting from a checkpoint,
540 * and in that case we should not do any modifications to the state.
542 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
544 // When restarting from a checkpoint, it can be appropriate to
545 // initialize ekind from quantities in the checkpoint. Otherwise,
546 // compute_globals must initialize ekind before the simulation
547 // starts/restarts. However, only the master rank knows what was
548 // found in the checkpoint file, so we have to communicate in
549 // order to coordinate the restart.
551 // TODO Consider removing this communication if/when checkpoint
552 // reading directly follows .tpr reading, because all ranks can
553 // agree on hasReadEkinState at that time.
554 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
557 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr->mpi_comm_mygroup);
559 if (hasReadEkinState)
561 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
564 unsigned int cglo_flags =
565 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
566 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
568 bSumEkinhOld = FALSE;
570 t_vcm vcm(top_global->groups, *ir);
571 reportComRemovalInfo(fplog, vcm);
573 /* To minimize communication, compute_globals computes the COM velocity
574 * and the kinetic energy for the velocities without COM motion removed.
575 * Thus to get the kinetic energy without the COM contribution, we need
576 * to call compute_globals twice.
578 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
580 unsigned int cglo_flags_iteration = cglo_flags;
581 if (bStopCM && cgloIteration == 0)
583 cglo_flags_iteration |= CGLO_STOPCM;
584 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
586 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
587 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, nullptr,
588 enerd, force_vir, shake_vir, total_vir, pres, constr, &nullSignaller,
589 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
591 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
593 if (cglo_flags_iteration & CGLO_STOPCM)
595 /* At initialization, do not pass x with acceleration-correction mode
596 * to avoid (incorrect) correction of the initial coordinates.
598 auto x = (vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION) ? ArrayRef<RVec>()
599 : makeArrayRef(state->x);
600 process_and_stopcm_grp(fplog, &vcm, *mdatoms, x, makeArrayRef(state->v));
601 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
604 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
605 makeConstArrayRef(state->x), state->box,
606 &shouldCheckNumberOfBondedInteractions);
607 if (ir->eI == eiVVAK)
609 /* a second call to get the half step temperature initialized as well */
610 /* we do the same call as above, but turn the pressure off -- internally to
611 compute_globals, this is recognized as a velocity verlet half-step
612 kinetic energy calculation. This minimized excess variables, but
613 perhaps loses some logic?*/
615 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
616 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, nullptr,
617 enerd, force_vir, shake_vir, total_vir, pres, constr, &nullSignaller,
618 state->box, nullptr, &bSumEkinhOld, cglo_flags & ~CGLO_PRESSURE);
621 /* Calculate the initial half step temperature, and save the ekinh_old */
622 if (startingBehavior == StartingBehavior::NewSimulation)
624 for (i = 0; (i < ir->opts.ngtc); i++)
626 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
630 /* need to make an initiation call to get the Trotter variables set, as well as other constants
631 for non-trotter temperature control */
632 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
636 if (!ir->bContinuation)
638 if (constr && ir->eConstrAlg == econtLINCS)
640 fprintf(fplog, "RMS relative constraint deviation after constraining: %.2e\n",
643 if (EI_STATE_VELOCITY(ir->eI))
645 real temp = enerd->term[F_TEMP];
648 /* Result of Ekin averaged over velocities of -half
649 * and +half step, while we only have -half step here.
653 fprintf(fplog, "Initial temperature: %g K\n", temp);
658 fprintf(stderr, "starting mdrun '%s'\n", *(top_global->name));
661 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
665 sprintf(tbuf, "%s", "infinite");
667 if (ir->init_step > 0)
669 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
670 gmx_step_str(ir->init_step + ir->nsteps, sbuf), tbuf,
671 gmx_step_str(ir->init_step, sbuf2), ir->init_step * ir->delta_t);
675 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
677 fprintf(fplog, "\n");
680 walltime_accounting_start_time(walltime_accounting);
681 wallcycle_start(wcycle, ewcRUN);
682 print_start(fplog, cr, walltime_accounting, "mdrun");
685 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
686 int chkpt_ret = fcCheckPointParallel(cr->nodeid, NULL, 0);
689 gmx_fatal(3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0);
693 /***********************************************************
697 ************************************************************/
700 /* Skip the first Nose-Hoover integration when we get the state from tpx */
701 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
702 bSumEkinhOld = FALSE;
704 bNeedRepartition = FALSE;
706 step = ir->init_step;
709 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
710 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
711 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
712 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
714 auto checkpointHandler = std::make_unique<CheckpointHandler>(
715 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]), simulationsShareState,
716 ir->nstlist == 0, MASTER(cr), mdrunOptions.writeConfout,
717 mdrunOptions.checkpointOptions.period);
719 const bool resetCountersIsLocal = true;
720 auto resetHandler = std::make_unique<ResetHandler>(
721 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
722 !resetCountersIsLocal, ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
723 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
725 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
727 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
729 logInitialMultisimStatus(ms, cr, mdlog, simulationsShareState, ir->nsteps, ir->init_step);
732 /* and stop now if we should */
733 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
737 /* Determine if this is a neighbor search step */
738 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
740 if (bPMETune && bNStList)
742 // This has to be here because PME load balancing is called so early.
743 // TODO: Move to after all booleans are defined.
744 if (useGpuForUpdate && !bFirstStep)
746 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
747 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
749 /* PME grid + cut-off optimization with GPUs or PME nodes */
750 pme_loadbal_do(pme_loadbal, cr, (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
751 fplog, mdlog, *ir, fr, state->box, state->x, wcycle, step, step_rel,
752 &bPMETunePrinting, simulationWork.useGpuPmePpCommunication);
755 wallcycle_start(wcycle, ewcSTEP);
757 bLastStep = (step_rel == ir->nsteps);
758 t = t0 + step * ir->delta_t;
760 // TODO Refactor this, so that nstfep does not need a default value of zero
761 if (ir->efep != efepNO || ir->bSimTemp)
763 /* find and set the current lambdas */
764 state->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
766 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
767 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
768 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
772 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
773 && do_per_step(step, replExParams.exchangeInterval));
775 if (doSimulatedAnnealing)
777 update_annealing_target_temp(ir, t, &upd);
780 /* Stop Center of Mass motion */
781 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
783 /* Determine whether or not to do Neighbour Searching */
784 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
786 /* Note that the stopHandler will cause termination at nstglobalcomm
787 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
788 * nstpcouple steps, we have computed the half-step kinetic energy
789 * of the previous step and can always output energies at the last step.
791 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
793 /* do_log triggers energy and virial calculation. Because this leads
794 * to different code paths, forces can be different. Thus for exact
795 * continuation we should avoid extra log output.
796 * Note that the || bLastStep can result in non-exact continuation
797 * beyond the last step. But we don't consider that to be an issue.
799 do_log = (do_per_step(step, ir->nstlog)
800 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
801 do_verbose = mdrunOptions.verbose
802 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
804 if (useGpuForUpdate && !bFirstStep && bNS)
806 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
807 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
808 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
809 // Copy coordinate from the GPU when needed at the search step.
810 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
811 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
812 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
813 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
816 if (bNS && !(bFirstStep && ir->bContinuation))
818 bMasterState = FALSE;
819 /* Correct the new box if it is too skewed */
820 if (inputrecDynamicBox(ir))
822 if (correct_box(fplog, step, state->box))
825 // If update is offloaded, it should be informed about the box size change
828 integrator->setPbc(PbcType::Xyz, state->box);
832 if (DOMAINDECOMP(cr) && bMasterState)
834 dd_collect_state(cr->dd, state, state_global);
837 if (DOMAINDECOMP(cr))
839 /* Repartition the domain decomposition */
840 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
841 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
842 fr, vsite, constr, nrnb, wcycle, do_verbose && !bPMETunePrinting);
843 shouldCheckNumberOfBondedInteractions = true;
844 upd.setNumAtoms(state->natoms);
848 // Allocate or re-size GPU halo exchange object, if necessary
849 if (bNS && havePPDomainDecomposition(cr) && simulationWork.useGpuHaloExchange
850 && useGpuForNonbonded && is1D(*cr->dd))
852 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
853 "GPU device manager has to be initialized to use GPU "
854 "version of halo exchange.");
855 // TODO remove need to pass local stream into GPU halo exchange - Issue #3093
856 constructGpuHaloExchange(mdlog, *cr, *fr->deviceStreamManager, wcycle);
859 if (MASTER(cr) && do_log)
861 gmx::EnergyOutput::printHeader(fplog, step,
862 t); /* can we improve the information printed here? */
865 if (ir->efep != efepNO)
867 update_mdatoms(mdatoms, state->lambda[efptMASS]);
873 /* We need the kinetic energy at minus the half step for determining
874 * the full step kinetic energy and possibly for T-coupling.*/
875 /* This may not be quite working correctly yet . . . . */
876 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
877 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, wcycle,
878 enerd, nullptr, nullptr, nullptr, nullptr, constr, &nullSignaller,
879 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
880 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
881 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
882 &top, makeConstArrayRef(state->x), state->box,
883 &shouldCheckNumberOfBondedInteractions);
885 clear_mat(force_vir);
887 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
889 /* Determine the energy and pressure:
890 * at nstcalcenergy steps and at energy output steps (set below).
892 if (EI_VV(ir->eI) && (!bInitStep))
894 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
895 bCalcVir = bCalcEnerStep
897 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
901 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
902 bCalcVir = bCalcEnerStep || (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
904 bCalcEner = bCalcEnerStep;
906 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
908 if (do_ene || do_log || bDoReplEx)
914 /* Do we need global communication ? */
915 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
916 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
918 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
919 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
920 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
924 /* Now is the time to relax the shells */
925 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
926 imdSession, pull_work, bNS, force_flags, &top, constr, enerd,
927 state->natoms, state->x.arrayRefWithPadding(),
928 state->v.arrayRefWithPadding(), state->box, state->lambda,
929 &state->hist, &f.view(), force_vir, mdatoms, nrnb, wcycle, shellfc,
930 fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
934 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
935 is updated (or the AWH update will be performed twice for one step when continuing).
936 It would be best to call this update function from do_md_trajectory_writing but that
937 would occur after do_force. One would have to divide the update_awh function into one
938 function applying the AWH force and one doing the AWH bias update. The update AWH
939 bias function could then be called after do_md_trajectory_writing (then containing
940 update_awh_history). The checkpointing will in the future probably moved to the start
941 of the md loop which will rid of this issue. */
942 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
944 awh->updateHistory(state_global->awhHistory.get());
947 /* The coordinates (x) are shifted (to get whole molecules)
949 * This is parallellized as well, and does communication too.
950 * Check comments in sim_util.c
952 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation, imdSession, pull_work, step,
953 nrnb, wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
954 &f.view(), force_vir, mdatoms, enerd, state->lambda, fr, runScheduleWork,
955 vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
956 (bNS ? GMX_FORCE_NS : 0) | force_flags, ddBalanceRegionHandler);
959 // VV integrators do not need the following velocity half step
960 // if it is the first step after starting from a checkpoint.
961 // That is, the half step is needed on all other steps, and
962 // also the first step when starting from a .tpr file.
963 if (EI_VV(ir->eI) && (!bFirstStep || startingBehavior == StartingBehavior::NewSimulation))
964 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
966 rvec* vbuf = nullptr;
968 wallcycle_start(wcycle, ewcUPDATE);
969 if (ir->eI == eiVV && bInitStep)
971 /* if using velocity verlet with full time step Ekin,
972 * take the first half step only to compute the
973 * virial for the first step. From there,
974 * revert back to the initial coordinates
975 * so that the input is actually the initial step.
977 snew(vbuf, state->natoms);
978 copy_rvecn(state->v.rvec_array(), vbuf, 0,
979 state->natoms); /* should make this better for parallelizing? */
983 /* this is for NHC in the Ekin(t+dt/2) version of vv */
984 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
985 trotter_seq, ettTSEQ1);
988 upd.update_coords(*ir, step, mdatoms, state, f.view().forceWithPadding(), fcdata, ekind,
989 M, etrtVELOCITY1, cr, constr != nullptr);
991 wallcycle_stop(wcycle, ewcUPDATE);
992 constrain_velocities(constr, do_log, do_ene, step, state, nullptr, bCalcVir, shake_vir);
993 wallcycle_start(wcycle, ewcUPDATE);
994 /* if VV, compute the pressure and constraints */
995 /* For VV2, we strictly only need this if using pressure
996 * control, but we really would like to have accurate pressures
998 * Think about ways around this in the future?
999 * For now, keep this choice in comments.
1001 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1002 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1004 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1005 if (bCalcEner && ir->eI == eiVVAK)
1007 bSumEkinhOld = TRUE;
1009 /* for vv, the first half of the integration actually corresponds to the previous step.
1010 So we need information from the last step in the first half of the integration */
1011 if (bGStat || do_per_step(step - 1, nstglobalcomm))
1013 wallcycle_stop(wcycle, ewcUPDATE);
1014 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1015 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, wcycle,
1016 enerd, force_vir, shake_vir, total_vir, pres, constr, &nullSignaller,
1017 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1018 (bGStat ? CGLO_GSTAT : 0) | (bCalcEner ? CGLO_ENERGY : 0)
1019 | (bTemp ? CGLO_TEMPERATURE : 0) | (bPres ? CGLO_PRESSURE : 0)
1020 | (bPres ? CGLO_CONSTRAINT : 0) | (bStopCM ? CGLO_STOPCM : 0)
1021 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1024 /* explanation of above:
1025 a) We compute Ekin at the full time step
1026 if 1) we are using the AveVel Ekin, and it's not the
1027 initial step, or 2) if we are using AveEkin, but need the full
1028 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1029 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1030 EkinAveVel because it's needed for the pressure */
1031 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1032 top_global, &top, makeConstArrayRef(state->x),
1033 state->box, &shouldCheckNumberOfBondedInteractions);
1036 process_and_stopcm_grp(fplog, &vcm, *mdatoms, makeArrayRef(state->x),
1037 makeArrayRef(state->v));
1038 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1040 wallcycle_start(wcycle, ewcUPDATE);
1042 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1047 m_add(force_vir, shake_vir,
1048 total_vir); /* we need the un-dispersion corrected total vir here */
1049 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1050 trotter_seq, ettTSEQ2);
1052 /* TODO This is only needed when we're about to write
1053 * a checkpoint, because we use it after the restart
1054 * (in a kludge?). But what should we be doing if
1055 * the startingBehavior is NewSimulation or bInitStep are true? */
1056 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1058 copy_mat(shake_vir, state->svir_prev);
1059 copy_mat(force_vir, state->fvir_prev);
1061 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1063 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1064 enerd->term[F_TEMP] =
1065 sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1066 enerd->term[F_EKIN] = trace(ekind->ekin);
1069 else if (bExchanged)
1071 wallcycle_stop(wcycle, ewcUPDATE);
1072 /* We need the kinetic energy at minus the half step for determining
1073 * the full step kinetic energy and possibly for T-coupling.*/
1074 /* This may not be quite working correctly yet . . . . */
1075 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1076 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, wcycle,
1077 enerd, nullptr, nullptr, nullptr, nullptr, constr, &nullSignaller,
1078 state->box, nullptr, &bSumEkinhOld, CGLO_GSTAT | CGLO_TEMPERATURE);
1079 wallcycle_start(wcycle, ewcUPDATE);
1082 /* if it's the initial step, we performed this first step just to get the constraint virial */
1083 if (ir->eI == eiVV && bInitStep)
1085 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1088 wallcycle_stop(wcycle, ewcUPDATE);
1091 /* compute the conserved quantity */
1094 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1097 last_ekin = enerd->term[F_EKIN];
1099 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1101 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1103 /* sum up the foreign kinetic energy and dK/dl terms for vv. currently done every step so that dhdl is correct in the .edr */
1104 if (ir->efep != efepNO)
1106 accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
1110 /* ######## END FIRST UPDATE STEP ############## */
1111 /* ######## If doing VV, we now have v(dt) ###### */
1114 /* perform extended ensemble sampling in lambda - we don't
1115 actually move to the new state before outputting
1116 statistics, but if performing simulated tempering, we
1117 do update the velocities and the tau_t. */
1119 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state,
1120 state->dfhist, step, state->v.rvec_array(), mdatoms);
1121 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1124 copy_df_history(state_global->dfhist, state->dfhist);
1128 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1129 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1130 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1131 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1132 || checkpointHandler->isCheckpointingStep()))
1134 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1135 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1137 // Copy velocities if needed for the output/checkpointing.
1138 // NOTE: Copy on the search steps is done at the beginning of the step.
1139 if (useGpuForUpdate && !bNS
1140 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1142 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1143 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1145 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1146 // and update is offloaded hence forces are kept on the GPU for update and have not been
1147 // already transferred in do_force().
1148 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1149 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1150 // prior to GPU update.
1151 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1152 // copy call in do_force(...).
1153 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1154 // on host after the D2H copy in do_force(...).
1155 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1156 && do_per_step(step, ir->nstfout))
1158 stateGpu->copyForcesFromGpu(f.view().force(), AtomLocality::Local);
1159 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1161 /* Now we have the energies and forces corresponding to the
1162 * coordinates at time t. We must output all of this before
1165 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state, state_global,
1166 observablesHistory, top_global, fr, outf, energyOutput, ekind,
1167 f.view().force(), checkpointHandler->isCheckpointingStep(),
1168 bRerunMD, bLastStep, mdrunOptions.writeConfout, bSumEkinhOld);
1169 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1170 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1172 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1173 if (startingBehavior != StartingBehavior::NewSimulation && bFirstStep
1174 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1176 copy_mat(state->svir_prev, shake_vir);
1177 copy_mat(state->fvir_prev, force_vir);
1180 stopHandler->setSignal();
1181 resetHandler->setSignal(walltime_accounting);
1183 if (bGStat || !PAR(cr))
1185 /* In parallel we only have to check for checkpointing in steps
1186 * where we do global communication,
1187 * otherwise the other nodes don't know.
1189 checkpointHandler->setSignal(walltime_accounting);
1192 /* ######### START SECOND UPDATE STEP ################# */
1194 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1195 controlled in preprocessing */
1197 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1199 gmx_bool bIfRandomize;
1200 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1201 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1202 if (constr && bIfRandomize)
1204 constrain_velocities(constr, do_log, do_ene, step, state, nullptr, false, nullptr);
1207 /* Box is changed in update() when we do pressure coupling,
1208 * but we should still use the old box for energy corrections and when
1209 * writing it to the energy file, so it matches the trajectory files for
1210 * the same timestep above. Make a copy in a separate array.
1212 copy_mat(state->box, lastbox);
1216 wallcycle_start(wcycle, ewcUPDATE);
1217 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1220 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1221 /* We can only do Berendsen coupling after we have summed
1222 * the kinetic energy or virial. Since the happens
1223 * in global_state after update, we should only do it at
1224 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1229 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1230 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1235 /* velocity half-step update */
1236 upd.update_coords(*ir, step, mdatoms, state, f.view().forceWithPadding(), fcdata, ekind,
1237 M, etrtVELOCITY2, cr, constr != nullptr);
1240 /* Above, initialize just copies ekinh into ekin,
1241 * it doesn't copy position (for VV),
1242 * and entire integrator for MD.
1245 if (ir->eI == eiVVAK)
1247 cbuf.resize(state->x.size());
1248 std::copy(state->x.begin(), state->x.end(), cbuf.begin());
1251 /* With leap-frog type integrators we compute the kinetic energy
1252 * at a whole time step as the average of the half-time step kinetic
1253 * energies of two subsequent steps. Therefore we need to compute the
1254 * half step kinetic energy also if we need energies at the next step.
1256 const bool needHalfStepKineticEnergy =
1257 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1259 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1260 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1261 const bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN
1262 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1264 if (useGpuForUpdate)
1266 if (bNS && (bFirstStep || DOMAINDECOMP(cr)))
1268 integrator->set(stateGpu->getCoordinates(), stateGpu->getVelocities(),
1269 stateGpu->getForces(), top.idef, *mdatoms, ekind->ngtc);
1271 // Copy data to the GPU after buffers might have being reinitialized
1272 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1273 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1276 // If the buffer ops were not offloaded this step, the forces are on the host and have to be copied
1277 if (!runScheduleWork->stepWork.useGpuFBufferOps)
1279 stateGpu->copyForcesToGpu(f.view().force(), AtomLocality::Local);
1282 const bool doTemperatureScaling =
1283 (ir->etc != etcNO && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1285 // This applies Leap-Frog, LINCS and SETTLE in succession
1286 integrator->integrate(stateGpu->getForcesReadyOnDeviceEvent(
1287 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1288 ir->delta_t, true, bCalcVir, shake_vir, doTemperatureScaling,
1289 ekind->tcstat, doParrinelloRahman, ir->nstpcouple * ir->delta_t, M);
1291 // Copy velocities D2H after update if:
1292 // - Globals are computed this step (includes the energy output steps).
1293 // - Temperature is needed for the next step.
1294 if (bGStat || needHalfStepKineticEnergy)
1296 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1297 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1302 upd.update_coords(*ir, step, mdatoms, state, f.view().forceWithPadding(), fcdata, ekind,
1303 M, etrtPOSITION, cr, constr != nullptr);
1305 wallcycle_stop(wcycle, ewcUPDATE);
1307 constrain_coordinates(constr, do_log, do_ene, step, state,
1308 upd.xp()->arrayRefWithPadding(), &dvdl_constr, bCalcVir, shake_vir);
1310 upd.update_sd_second_half(*ir, step, &dvdl_constr, mdatoms, state, cr, nrnb, wcycle,
1311 constr, do_log, do_ene);
1312 upd.finish_update(*ir, mdatoms, state, wcycle, constr != nullptr);
1315 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1317 updatePrevStepPullCom(pull_work, state);
1320 if (ir->eI == eiVVAK)
1322 /* erase F_EKIN and F_TEMP here? */
1323 /* just compute the kinetic energy at the half step to perform a trotter step */
1324 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1325 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, wcycle, enerd,
1326 force_vir, shake_vir, total_vir, pres, constr, &nullSignaller, lastbox,
1327 nullptr, &bSumEkinhOld, (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE);
1328 wallcycle_start(wcycle, ewcUPDATE);
1329 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1330 /* now we know the scaling, we can compute the positions again */
1331 std::copy(cbuf.begin(), cbuf.end(), state->x.begin());
1333 upd.update_coords(*ir, step, mdatoms, state, f.view().forceWithPadding(), fcdata, ekind,
1334 M, etrtPOSITION, cr, constr != nullptr);
1335 wallcycle_stop(wcycle, ewcUPDATE);
1337 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1338 /* are the small terms in the shake_vir here due
1339 * to numerical errors, or are they important
1340 * physically? I'm thinking they are just errors, but not completely sure.
1341 * For now, will call without actually constraining, constr=NULL*/
1342 upd.finish_update(*ir, mdatoms, state, wcycle, false);
1346 /* this factor or 2 correction is necessary
1347 because half of the constraint force is removed
1348 in the vv step, so we have to double it. See
1349 the Issue #1255. It is not yet clear
1350 if the factor of 2 is exact, or just a very
1351 good approximation, and this will be
1352 investigated. The next step is to see if this
1353 can be done adding a dhdl contribution from the
1354 rattle step, but this is somewhat more
1355 complicated with the current code. Will be
1356 investigated, hopefully for 4.6.3. However,
1357 this current solution is much better than
1358 having it completely wrong.
1360 enerd->term[F_DVDL_CONSTR] += 2 * dvdl_constr;
1364 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1367 if (vsite != nullptr)
1369 wallcycle_start(wcycle, ewcVSITECONSTR);
1370 vsite->construct(state->x, ir->delta_t, state->v, state->box);
1371 wallcycle_stop(wcycle, ewcVSITECONSTR);
1374 /* ############## IF NOT VV, Calculate globals HERE ############ */
1375 /* With Leap-Frog we can skip compute_globals at
1376 * non-communication steps, but we need to calculate
1377 * the kinetic energy one step before communication.
1380 // Organize to do inter-simulation signalling on steps if
1381 // and when algorithms require it.
1382 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1384 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1386 // Copy coordinates when needed to stop the CM motion.
1387 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1389 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1390 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1392 // Since we're already communicating at this step, we
1393 // can propagate intra-simulation signals. Note that
1394 // check_nstglobalcomm has the responsibility for
1395 // choosing the value of nstglobalcomm that is one way
1396 // bGStat becomes true, so we can't get into a
1397 // situation where e.g. checkpointing can't be
1399 bool doIntraSimSignal = true;
1400 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1402 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1403 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm,
1404 wcycle, enerd, force_vir, shake_vir, total_vir, pres, constr,
1405 &signaller, lastbox, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1406 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1407 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1408 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1409 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1410 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1412 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1413 top_global, &top, makeConstArrayRef(state->x),
1414 state->box, &shouldCheckNumberOfBondedInteractions);
1415 if (!EI_VV(ir->eI) && bStopCM)
1417 process_and_stopcm_grp(fplog, &vcm, *mdatoms, makeArrayRef(state->x),
1418 makeArrayRef(state->v));
1419 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1421 // TODO: The special case of removing CM motion should be dealt more gracefully
1422 if (useGpuForUpdate)
1424 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1425 // Here we block until the H2D copy completes because event sync with the
1426 // force kernels that use the coordinates on the next steps is not implemented
1427 // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1428 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1429 // If the COM removal changed the velocities on the CPU, this has to be accounted for.
1430 if (vcm.mode != ecmNO)
1432 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1439 /* ############# END CALC EKIN AND PRESSURE ################# */
1441 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1442 the virial that should probably be addressed eventually. state->veta has better properies,
1443 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1444 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1446 if (ir->efep != efepNO && !EI_VV(ir->eI))
1448 /* Sum up the foreign energy and dK/dl terms for md and sd.
1449 Currently done every step so that dH/dl is correct in the .edr */
1450 accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
1453 update_pcouple_after_coordinates(fplog, step, ir, mdatoms, pres, force_vir, shake_vir,
1454 pressureCouplingMu, state, nrnb, upd.deform(), !useGpuForUpdate);
1456 const bool doBerendsenPressureCoupling =
1457 (inputrec->epc == epcBERENDSEN && do_per_step(step, inputrec->nstpcouple));
1458 if (useGpuForUpdate && (doBerendsenPressureCoupling || doParrinelloRahman))
1460 integrator->scaleCoordinates(pressureCouplingMu);
1461 integrator->setPbc(PbcType::Xyz, state->box);
1464 /* ################# END UPDATE STEP 2 ################# */
1465 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1467 /* The coordinates (x) were unshifted in update */
1470 /* We will not sum ekinh_old,
1471 * so signal that we still have to do it.
1473 bSumEkinhOld = TRUE;
1478 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1480 /* use the directly determined last velocity, not actually the averaged half steps */
1481 if (bTrotter && ir->eI == eiVV)
1483 enerd->term[F_EKIN] = last_ekin;
1485 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1487 if (integratorHasConservedEnergyQuantity(ir))
1491 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1495 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1498 /* ######### END PREPARING EDR OUTPUT ########### */
1504 if (fplog && do_log && bDoExpanded)
1506 /* only needed if doing expanded ensemble */
1507 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals,
1508 ir->bSimTemp ? ir->simtempvals : nullptr,
1509 state_global->dfhist, state->fep_state, ir->nstlog, step);
1513 energyOutput.addDataAtEnergyStep(
1514 bDoDHDL, bCalcEnerStep, t, mdatoms->tmass, enerd, ir->fepvals,
1515 ir->expandedvals, lastbox,
1516 PTCouplingArrays{ state->boxv, state->nosehoover_xi, state->nosehoover_vxi,
1517 state->nhpres_xi, state->nhpres_vxi },
1518 state->fep_state, shake_vir, force_vir, total_vir, pres, ekind, mu_tot, constr);
1522 energyOutput.recordNonEnergyStep();
1525 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1526 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1528 if (doSimulatedAnnealing)
1530 gmx::EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups,
1533 if (do_log || do_ene || do_dr || do_or)
1535 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1536 do_log ? fplog : nullptr, step, t,
1537 fr->fcdata.get(), awh.get());
1539 if (do_log && ir->bDoAwh && awh->hasFepLambdaDimension())
1541 const bool isInitialOutput = false;
1542 printLambdaStateToLog(fplog, state->lambda, isInitialOutput);
1547 pull_print_output(pull_work, step, t);
1550 if (do_per_step(step, ir->nstlog))
1552 if (fflush(fplog) != 0)
1554 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1560 /* Have to do this part _after_ outputting the logfile and the edr file */
1561 /* Gets written into the state at the beginning of next loop*/
1562 state->fep_state = lamnew;
1564 else if (ir->bDoAwh && awh->needForeignEnergyDifferences(step))
1566 state->fep_state = awh->fepLambdaState();
1568 /* Print the remaining wall clock time for the run */
1569 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1573 fprintf(stderr, "\n");
1575 print_time(stderr, walltime_accounting, step, ir, cr);
1578 /* Ion/water position swapping.
1579 * Not done in last step since trajectory writing happens before this call
1580 * in the MD loop and exchanges would be lost anyway. */
1581 bNeedRepartition = FALSE;
1582 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep && do_per_step(step, ir->swap->nstswap))
1585 do_swapcoords(cr, step, t, ir, swap, wcycle, as_rvec_array(state->x.data()),
1586 state->box, MASTER(cr) && mdrunOptions.verbose, bRerunMD);
1588 if (bNeedRepartition && DOMAINDECOMP(cr))
1590 dd_collect_state(cr->dd, state, state_global);
1594 /* Replica exchange */
1598 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1601 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1603 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1, state_global, *top_global, ir,
1604 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
1605 nrnb, wcycle, FALSE);
1606 shouldCheckNumberOfBondedInteractions = true;
1607 upd.setNumAtoms(state->natoms);
1613 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1614 /* With all integrators, except VV, we need to retain the pressure
1615 * at the current step for coupling at the next step.
1617 if ((state->flags & (1U << estPRES_PREV))
1618 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1620 /* Store the pressure in t_state for pressure coupling
1621 * at the next MD step.
1623 copy_mat(pres, state->pres_prev);
1626 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1628 if ((membed != nullptr) && (!bLastStep))
1630 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1633 cycles = wallcycle_stop(wcycle, ewcSTEP);
1634 if (DOMAINDECOMP(cr) && wcycle)
1636 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1639 /* increase the MD step number */
1643 resetHandler->resetCounters(step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb,
1644 fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1646 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1647 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1649 /* End of main MD loop */
1651 /* Closing TNG files can include compressing data. Therefore it is good to do that
1652 * before stopping the time measurements. */
1653 mdoutf_tng_close(outf);
1655 /* Stop measuring walltime */
1656 walltime_accounting_end_time(walltime_accounting);
1658 if (!thisRankHasDuty(cr, DUTY_PME))
1660 /* Tell the PME only node to finish */
1661 gmx_pme_send_finish(cr);
1666 if (ir->nstcalcenergy > 0)
1668 gmx::EnergyOutput::printAnnealingTemperatures(fplog, groups, &(ir->opts));
1669 energyOutput.printAverages(fplog, groups);
1676 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1679 done_shellfc(fplog, shellfc, step_rel);
1681 if (useReplicaExchange && MASTER(cr))
1683 print_replica_exchange_statistics(fplog, repl_ex);
1686 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1688 global_stat_destroy(gstat);