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
54 #include "gromacs/awh/awh.h"
55 #include "gromacs/commandline/filenm.h"
56 #include "gromacs/domdec/collect.h"
57 #include "gromacs/domdec/dlbtiming.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_network.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/domdec/gpuhaloexchange.h"
62 #include "gromacs/domdec/mdsetup.h"
63 #include "gromacs/domdec/partition.h"
64 #include "gromacs/essentialdynamics/edsam.h"
65 #include "gromacs/ewald/pme_load_balancing.h"
66 #include "gromacs/ewald/pme_pp.h"
67 #include "gromacs/fileio/trxio.h"
68 #include "gromacs/gmxlib/network.h"
69 #include "gromacs/gmxlib/nrnb.h"
70 #include "gromacs/gpu_utils/gpu_utils.h"
71 #include "gromacs/imd/imd.h"
72 #include "gromacs/listed_forces/manage_threading.h"
73 #include "gromacs/math/functions.h"
74 #include "gromacs/math/utilities.h"
75 #include "gromacs/math/vec.h"
76 #include "gromacs/math/vectypes.h"
77 #include "gromacs/mdlib/checkpointhandler.h"
78 #include "gromacs/mdlib/compute_io.h"
79 #include "gromacs/mdlib/constr.h"
80 #include "gromacs/mdlib/ebin.h"
81 #include "gromacs/mdlib/enerdata_utils.h"
82 #include "gromacs/mdlib/energyoutput.h"
83 #include "gromacs/mdlib/expanded.h"
84 #include "gromacs/mdlib/force.h"
85 #include "gromacs/mdlib/force_flags.h"
86 #include "gromacs/mdlib/forcerec.h"
87 #include "gromacs/mdlib/md_support.h"
88 #include "gromacs/mdlib/mdatoms.h"
89 #include "gromacs/mdlib/mdoutf.h"
90 #include "gromacs/mdlib/membed.h"
91 #include "gromacs/mdlib/resethandler.h"
92 #include "gromacs/mdlib/sighandler.h"
93 #include "gromacs/mdlib/simulationsignal.h"
94 #include "gromacs/mdlib/stat.h"
95 #include "gromacs/mdlib/stophandler.h"
96 #include "gromacs/mdlib/tgroup.h"
97 #include "gromacs/mdlib/trajectory_writing.h"
98 #include "gromacs/mdlib/update.h"
99 #include "gromacs/mdlib/update_constrain_gpu.h"
100 #include "gromacs/mdlib/vcm.h"
101 #include "gromacs/mdlib/vsite.h"
102 #include "gromacs/mdrunutility/handlerestart.h"
103 #include "gromacs/mdrunutility/multisim.h"
104 #include "gromacs/mdrunutility/printtime.h"
105 #include "gromacs/mdtypes/awh_history.h"
106 #include "gromacs/mdtypes/awh_params.h"
107 #include "gromacs/mdtypes/commrec.h"
108 #include "gromacs/mdtypes/df_history.h"
109 #include "gromacs/mdtypes/energyhistory.h"
110 #include "gromacs/mdtypes/fcdata.h"
111 #include "gromacs/mdtypes/forcerec.h"
112 #include "gromacs/mdtypes/group.h"
113 #include "gromacs/mdtypes/inputrec.h"
114 #include "gromacs/mdtypes/interaction_const.h"
115 #include "gromacs/mdtypes/md_enums.h"
116 #include "gromacs/mdtypes/mdatom.h"
117 #include "gromacs/mdtypes/mdrunoptions.h"
118 #include "gromacs/mdtypes/observableshistory.h"
119 #include "gromacs/mdtypes/pullhistory.h"
120 #include "gromacs/mdtypes/simulation_workload.h"
121 #include "gromacs/mdtypes/state.h"
122 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
123 #include "gromacs/modularsimulator/energyelement.h"
124 #include "gromacs/nbnxm/gpu_data_mgmt.h"
125 #include "gromacs/nbnxm/nbnxm.h"
126 #include "gromacs/pbcutil/mshift.h"
127 #include "gromacs/pbcutil/pbc.h"
128 #include "gromacs/pulling/output.h"
129 #include "gromacs/pulling/pull.h"
130 #include "gromacs/swap/swapcoords.h"
131 #include "gromacs/timing/wallcycle.h"
132 #include "gromacs/timing/walltime_accounting.h"
133 #include "gromacs/topology/atoms.h"
134 #include "gromacs/topology/idef.h"
135 #include "gromacs/topology/mtop_util.h"
136 #include "gromacs/topology/topology.h"
137 #include "gromacs/trajectory/trajectoryframe.h"
138 #include "gromacs/utility/basedefinitions.h"
139 #include "gromacs/utility/cstringutil.h"
140 #include "gromacs/utility/fatalerror.h"
141 #include "gromacs/utility/logger.h"
142 #include "gromacs/utility/real.h"
143 #include "gromacs/utility/smalloc.h"
145 #include "legacysimulator.h"
146 #include "replicaexchange.h"
150 # include "corewrap.h"
153 using gmx::SimulationSignaller;
155 void gmx::LegacySimulator::do_md()
157 // TODO Historically, the EM and MD "integrators" used different
158 // names for the t_inputrec *parameter, but these must have the
159 // same name, now that it's a member of a struct. We use this ir
160 // alias to avoid a large ripple of nearly useless changes.
161 // t_inputrec is being replaced by IMdpOptionsProvider, so this
162 // will go away eventually.
163 t_inputrec* ir = inputrec;
164 int64_t step, step_rel;
165 double t, t0 = ir->init_t, lam0[efptNR];
166 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
167 gmx_bool bNS, bNStList, bStopCM, bFirstStep, bInitStep, bLastStep = FALSE;
168 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
169 gmx_bool do_ene, do_log, do_verbose;
170 gmx_bool bMasterState;
171 unsigned int force_flags;
172 tensor force_vir = { { 0 } }, shake_vir = { { 0 } }, total_vir = { { 0 } }, tmp_vir = { { 0 } },
176 matrix pressureCouplingMu, M;
177 gmx_repl_ex_t repl_ex = nullptr;
179 PaddedHostVector<gmx::RVec> f{};
180 gmx_global_stat_t gstat;
181 t_graph* graph = nullptr;
182 gmx_shellfc_t* shellfc;
183 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
184 gmx_bool bTemp, bPres, bTrotter;
186 std::vector<RVec> cbuf;
192 real saved_conserved_quantity = 0;
195 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
197 /* PME load balancing data for GPU kernels */
198 gmx_bool bPMETune = FALSE;
199 gmx_bool bPMETunePrinting = FALSE;
201 bool bInteractiveMDstep = false;
203 /* Domain decomposition could incorrectly miss a bonded
204 interaction, but checking for that requires a global
205 communication stage, which does not otherwise happen in DD
206 code. So we do that alongside the first global energy reduction
207 after a new DD is made. These variables handle whether the
208 check happens, and the result it returns. */
209 bool shouldCheckNumberOfBondedInteractions = false;
210 int totalNumberOfBondedInteractions = -1;
212 SimulationSignals signals;
213 // Most global communnication stages don't propagate mdrun
214 // signals, and will use this object to achieve that.
215 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
217 if (!mdrunOptions.writeConfout)
219 // This is on by default, and the main known use case for
220 // turning it off is for convenience in benchmarking, which is
221 // something that should not show up in the general user
226 "The -noconfout functionality is deprecated, and may be removed in a "
230 /* md-vv uses averaged full step velocities for T-control
231 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
232 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
233 bTrotter = (EI_VV(ir->eI)
234 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
236 const bool bRerunMD = false;
238 int nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
239 bGStatEveryStep = (nstglobalcomm == 1);
241 SimulationGroups* groups = &top_global->groups;
243 std::unique_ptr<EssentialDynamics> ed = nullptr;
244 if (opt2bSet("-ei", nfile, fnm))
246 /* Initialize essential dynamics sampling */
247 ed = init_edsam(mdlog, opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm), top_global,
248 ir, cr, constr, state_global, observablesHistory, oenv, startingBehavior);
250 else if (observablesHistory->edsamHistory)
253 "The checkpoint is from a run with essential dynamics sampling, "
254 "but the current run did not specify the -ei option. "
255 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
258 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
259 Update upd(ir, deform);
260 const bool doSimulatedAnnealing = initSimulatedAnnealing(ir, &upd);
261 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
263 bool simulationsShareState = false;
264 int nstSignalComm = nstglobalcomm;
266 // TODO This implementation of ensemble orientation restraints is nasty because
267 // a user can't just do multi-sim with single-sim orientation restraints.
268 bool usingEnsembleRestraints =
269 (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
270 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
272 // Replica exchange, ensemble restraints and AWH need all
273 // simulations to remain synchronized, so they need
274 // checkpoints and stop conditions to act on the same step, so
275 // the propagation of such signals must take place between
276 // simulations, not just within simulations.
277 // TODO: Make algorithm initializers set these flags.
278 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
280 if (simulationsShareState)
282 // Inter-simulation signal communication does not need to happen
283 // often, so we use a minimum of 200 steps to reduce overhead.
284 const int c_minimumInterSimulationSignallingInterval = 200;
285 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
290 if (startingBehavior != StartingBehavior::RestartWithAppending)
292 pleaseCiteCouplingAlgorithms(fplog, *ir);
295 init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, ir,
296 top_global, oenv, wcycle, startingBehavior, simulationsShareState, ms);
297 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
298 mdoutf_get_fp_dhdl(outf), false, startingBehavior, mdModulesNotifier);
300 gstat = global_stat_init(ir);
302 /* Check for polarizable models and flexible constraints */
303 shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
304 ir->nstcalcenergy, DOMAINDECOMP(cr));
307 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
308 if ((io > 2000) && MASTER(cr))
310 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
314 // Local state only becomes valid now.
315 std::unique_ptr<t_state> stateInstance;
319 auto mdatoms = mdAtoms->mdatoms();
321 std::unique_ptr<UpdateConstrainGpu> integrator;
323 if (DOMAINDECOMP(cr))
325 dd_init_local_top(*top_global, &top);
327 stateInstance = std::make_unique<t_state>();
328 state = stateInstance.get();
329 dd_init_local_state(cr->dd, state_global, state);
331 /* Distribute the charge groups over the nodes from the master node */
332 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
333 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
334 nrnb, nullptr, FALSE);
335 shouldCheckNumberOfBondedInteractions = true;
336 upd.setNumAtoms(state->natoms);
340 state_change_natoms(state_global, state_global->natoms);
341 f.resizeWithPadding(state_global->natoms);
342 /* Copy the pointer to the global state */
343 state = state_global;
345 /* Generate and initialize new topology */
346 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &graph, mdAtoms, constr, vsite, shellfc);
348 upd.setNumAtoms(state->natoms);
351 const auto& simulationWork = runScheduleWork->simulationWork;
352 const bool useGpuForPme = simulationWork.useGpuPme;
353 const bool useGpuForNonbonded = simulationWork.useGpuNonbonded;
354 const bool useGpuForBufferOps = simulationWork.useGpuBufferOps;
355 const bool useGpuForUpdate = simulationWork.useGpuUpdate;
357 StatePropagatorDataGpu* stateGpu = fr->stateGpu;
361 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
362 || constr->numConstraintsTotal() == 0,
363 "Constraints in domain decomposition are only supported with update "
364 "groups if using GPU update.\n");
365 GMX_RELEASE_ASSERT(ir->eConstrAlg != econtSHAKE || constr == nullptr
366 || constr->numConstraintsTotal() == 0,
367 "SHAKE is not supported with GPU update.");
368 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
369 "Either PME or short-ranged non-bonded interaction tasks must run on "
370 "the GPU to use GPU update.\n");
371 GMX_RELEASE_ASSERT(ir->eI == eiMD,
372 "Only the md integrator is supported with the GPU update.\n");
374 ir->etc != etcNOSEHOOVER,
375 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
376 GMX_RELEASE_ASSERT(ir->epc == epcNO || ir->epc == epcPARRINELLORAHMAN || ir->epc == epcBERENDSEN,
377 "Only Parrinello-Rahman and Berendsen pressure coupling are supported "
378 "with the GPU update.\n");
379 GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
380 "Virtual sites are not supported with the GPU update.\n");
381 GMX_RELEASE_ASSERT(ed == nullptr,
382 "Essential dynamics is not supported with the GPU update.\n");
383 GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(ir->pull),
384 "Constraints pulling is not supported with the GPU update.\n");
385 GMX_RELEASE_ASSERT(fcd->orires.nr == 0,
386 "Orientation restraints are not supported with the GPU update.\n");
389 || (!haveFreeEnergyType(*ir, efptBONDED) && !haveFreeEnergyType(*ir, efptMASS)),
390 "Free energy perturbation of masses and constraints are not supported with the GPU "
392 GMX_RELEASE_ASSERT(graph == nullptr, "The graph is not supported with GPU update.");
394 if (constr != nullptr && constr->numConstraintsTotal() > 0)
398 .appendText("Updating coordinates and applying constraints on the GPU.");
402 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
404 integrator = std::make_unique<UpdateConstrainGpu>(
405 *ir, *top_global, stateGpu->getUpdateStream(), stateGpu->xUpdatedOnDevice());
407 integrator->setPbc(PbcType::Xyz, state->box);
410 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
412 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
414 if ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
416 changePinningPolicy(&f, PinningPolicy::PinnedIfSupported);
420 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
423 // NOTE: The global state is no longer used at this point.
424 // But state_global is still used as temporary storage space for writing
425 // the global state to file and potentially for replica exchange.
426 // (Global topology should persist.)
428 update_mdatoms(mdatoms, state->lambda[efptMASS]);
432 /* Check nstexpanded here, because the grompp check was broken */
433 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
436 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
438 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
443 EnergyElement::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
446 preparePrevStepPullCom(ir, pull_work, mdatoms, state, state_global, cr,
447 startingBehavior != StartingBehavior::NewSimulation);
449 // TODO: Remove this by converting AWH into a ForceProvider
450 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms,
451 startingBehavior != StartingBehavior::NewSimulation,
452 shellfc != nullptr, opt2fn("-awh", nfile, fnm), pull_work);
454 if (useReplicaExchange && MASTER(cr))
456 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir, replExParams);
458 /* PME tuning is only supported in the Verlet scheme, with PME for
459 * Coulomb. It is not supported with only LJ PME. */
460 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
461 && ir->cutoff_scheme != ecutsGROUP);
463 pme_load_balancing_t* pme_loadbal = nullptr;
466 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata,
470 if (!ir->bContinuation)
472 if (state->flags & (1U << estV))
474 auto v = makeArrayRef(state->v);
475 /* Set the velocities of vsites, shells and frozen atoms to zero */
476 for (i = 0; i < mdatoms->homenr; i++)
478 if (mdatoms->ptype[i] == eptVSite || mdatoms->ptype[i] == eptShell)
482 else if (mdatoms->cFREEZE)
484 for (m = 0; m < DIM; m++)
486 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
497 /* Constrain the initial coordinates and velocities */
498 do_constrain_first(fplog, constr, ir, mdatoms, state->natoms, state->x.arrayRefWithPadding(),
499 state->v.arrayRefWithPadding(), state->box, state->lambda[efptBONDED]);
503 /* Construct the virtual sites for the initial configuration */
504 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, nullptr, top.idef.iparams,
505 top.idef.il, fr->pbcType, fr->bMolPBC, cr, state->box);
509 if (ir->efep != efepNO)
511 /* Set free energy calculation frequency as the greatest common
512 * denominator of nstdhdl and repl_ex_nst. */
513 nstfep = ir->fepvals->nstdhdl;
516 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
518 if (useReplicaExchange)
520 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
524 /* Be REALLY careful about what flags you set here. You CANNOT assume
525 * this is the first step, since we might be restarting from a checkpoint,
526 * and in that case we should not do any modifications to the state.
528 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
530 // When restarting from a checkpoint, it can be appropriate to
531 // initialize ekind from quantities in the checkpoint. Otherwise,
532 // compute_globals must initialize ekind before the simulation
533 // starts/restarts. However, only the master rank knows what was
534 // found in the checkpoint file, so we have to communicate in
535 // order to coordinate the restart.
537 // TODO Consider removing this communication if/when checkpoint
538 // reading directly follows .tpr reading, because all ranks can
539 // agree on hasReadEkinState at that time.
540 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
543 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr);
545 if (hasReadEkinState)
547 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
550 unsigned int cglo_flags =
551 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
552 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
554 bSumEkinhOld = FALSE;
556 t_vcm vcm(top_global->groups, *ir);
557 reportComRemovalInfo(fplog, vcm);
559 /* To minimize communication, compute_globals computes the COM velocity
560 * and the kinetic energy for the velocities without COM motion removed.
561 * Thus to get the kinetic energy without the COM contribution, we need
562 * to call compute_globals twice.
564 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
566 unsigned int cglo_flags_iteration = cglo_flags;
567 if (bStopCM && cgloIteration == 0)
569 cglo_flags_iteration |= CGLO_STOPCM;
570 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
572 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
573 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
574 nrnb, &vcm, nullptr, enerd, force_vir, shake_vir, total_vir, pres, constr,
575 &nullSignaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
577 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
579 if (cglo_flags_iteration & CGLO_STOPCM)
581 /* At initialization, do not pass x with acceleration-correction mode
582 * to avoid (incorrect) correction of the initial coordinates.
584 auto x = (vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION) ? ArrayRef<RVec>()
585 : makeArrayRef(state->x);
586 process_and_stopcm_grp(fplog, &vcm, *mdatoms, x, makeArrayRef(state->v));
587 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
590 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
591 makeConstArrayRef(state->x), state->box,
592 &shouldCheckNumberOfBondedInteractions);
593 if (ir->eI == eiVVAK)
595 /* a second call to get the half step temperature initialized as well */
596 /* we do the same call as above, but turn the pressure off -- internally to
597 compute_globals, this is recognized as a velocity verlet half-step
598 kinetic energy calculation. This minimized excess variables, but
599 perhaps loses some logic?*/
601 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
602 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
603 nrnb, &vcm, nullptr, enerd, force_vir, shake_vir, total_vir, pres, constr,
604 &nullSignaller, state->box, nullptr, &bSumEkinhOld, cglo_flags & ~CGLO_PRESSURE);
607 /* Calculate the initial half step temperature, and save the ekinh_old */
608 if (startingBehavior == StartingBehavior::NewSimulation)
610 for (i = 0; (i < ir->opts.ngtc); i++)
612 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
616 /* need to make an initiation call to get the Trotter variables set, as well as other constants
617 for non-trotter temperature control */
618 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
622 if (!ir->bContinuation)
624 if (constr && ir->eConstrAlg == econtLINCS)
626 fprintf(fplog, "RMS relative constraint deviation after constraining: %.2e\n",
629 if (EI_STATE_VELOCITY(ir->eI))
631 real temp = enerd->term[F_TEMP];
634 /* Result of Ekin averaged over velocities of -half
635 * and +half step, while we only have -half step here.
639 fprintf(fplog, "Initial temperature: %g K\n", temp);
644 fprintf(stderr, "starting mdrun '%s'\n", *(top_global->name));
647 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
651 sprintf(tbuf, "%s", "infinite");
653 if (ir->init_step > 0)
655 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
656 gmx_step_str(ir->init_step + ir->nsteps, sbuf), tbuf,
657 gmx_step_str(ir->init_step, sbuf2), ir->init_step * ir->delta_t);
661 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
663 fprintf(fplog, "\n");
666 walltime_accounting_start_time(walltime_accounting);
667 wallcycle_start(wcycle, ewcRUN);
668 print_start(fplog, cr, walltime_accounting, "mdrun");
671 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
672 int chkpt_ret = fcCheckPointParallel(cr->nodeid, NULL, 0);
675 gmx_fatal(3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0);
679 /***********************************************************
683 ************************************************************/
686 /* Skip the first Nose-Hoover integration when we get the state from tpx */
687 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
688 bSumEkinhOld = FALSE;
690 bNeedRepartition = FALSE;
692 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
693 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
694 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
695 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
697 auto checkpointHandler = std::make_unique<CheckpointHandler>(
698 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]), simulationsShareState,
699 ir->nstlist == 0, MASTER(cr), mdrunOptions.writeConfout,
700 mdrunOptions.checkpointOptions.period);
702 const bool resetCountersIsLocal = true;
703 auto resetHandler = std::make_unique<ResetHandler>(
704 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
705 !resetCountersIsLocal, ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
706 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
708 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
710 step = ir->init_step;
713 // TODO extract this to new multi-simulation module
714 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
716 if (!multisim_int_all_are_equal(ms, ir->nsteps))
718 GMX_LOG(mdlog.warning)
720 "Note: The number of steps is not consistent across multi "
722 "but we are proceeding anyway!");
724 if (!multisim_int_all_are_equal(ms, ir->init_step))
726 if (simulationsShareState)
731 "The initial step is not consistent across multi simulations which "
738 GMX_LOG(mdlog.warning)
740 "Note: The initial step is not consistent across multi "
742 "but we are proceeding anyway!");
747 /* and stop now if we should */
748 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
752 /* Determine if this is a neighbor search step */
753 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
755 if (bPMETune && bNStList)
757 // This has to be here because PME load balancing is called so early.
758 // TODO: Move to after all booleans are defined.
759 if (useGpuForUpdate && !bFirstStep)
761 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
762 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
764 /* PME grid + cut-off optimization with GPUs or PME nodes */
765 pme_loadbal_do(pme_loadbal, cr, (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
766 fplog, mdlog, *ir, fr, state->box, state->x, wcycle, step, step_rel,
767 &bPMETunePrinting, simulationWork.useGpuPmePpCommunication);
770 wallcycle_start(wcycle, ewcSTEP);
772 bLastStep = (step_rel == ir->nsteps);
773 t = t0 + step * ir->delta_t;
775 // TODO Refactor this, so that nstfep does not need a default value of zero
776 if (ir->efep != efepNO || ir->bSimTemp)
778 /* find and set the current lambdas */
779 setCurrentLambdasLocal(step, ir->fepvals, lam0, state->lambda, state->fep_state);
781 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
782 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
783 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
784 && (step > 0) && (startingBehavior == StartingBehavior::NewSimulation));
787 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
788 && do_per_step(step, replExParams.exchangeInterval));
790 if (doSimulatedAnnealing)
792 update_annealing_target_temp(ir, t, &upd);
795 /* Stop Center of Mass motion */
796 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
798 /* Determine whether or not to do Neighbour Searching */
799 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
801 /* Note that the stopHandler will cause termination at nstglobalcomm
802 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
803 * nstpcouple steps, we have computed the half-step kinetic energy
804 * of the previous step and can always output energies at the last step.
806 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
808 /* do_log triggers energy and virial calculation. Because this leads
809 * to different code paths, forces can be different. Thus for exact
810 * continuation we should avoid extra log output.
811 * Note that the || bLastStep can result in non-exact continuation
812 * beyond the last step. But we don't consider that to be an issue.
814 do_log = (do_per_step(step, ir->nstlog)
815 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
816 do_verbose = mdrunOptions.verbose
817 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
819 if (useGpuForUpdate && !bFirstStep && bNS)
821 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
822 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
823 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
824 // Copy coordinate from the GPU when needed at the search step.
825 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
826 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
827 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
828 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
831 if (bNS && !(bFirstStep && ir->bContinuation))
833 bMasterState = FALSE;
834 /* Correct the new box if it is too skewed */
835 if (inputrecDynamicBox(ir))
837 if (correct_box(fplog, step, state->box, graph))
840 // If update is offloaded, it should be informed about the box size change
843 integrator->setPbc(PbcType::Xyz, state->box);
847 if (DOMAINDECOMP(cr) && bMasterState)
849 dd_collect_state(cr->dd, state, state_global);
852 if (DOMAINDECOMP(cr))
854 /* Repartition the domain decomposition */
855 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
856 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
857 fr, vsite, constr, nrnb, wcycle, do_verbose && !bPMETunePrinting);
858 shouldCheckNumberOfBondedInteractions = true;
859 upd.setNumAtoms(state->natoms);
861 // Allocate or re-size GPU halo exchange object, if necessary
862 if (havePPDomainDecomposition(cr) && simulationWork.useGpuHaloExchange
863 && useGpuForNonbonded && is1D(*cr->dd))
865 // TODO remove need to pass local stream into GPU halo exchange - Redmine #3093
867 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local);
868 void* streamNonLocal = Nbnxm::gpu_get_command_stream(
869 fr->nbv->gpu_nbv, InteractionLocality::NonLocal);
870 constructGpuHaloExchange(mdlog, *cr, streamLocal, streamNonLocal);
875 if (MASTER(cr) && do_log)
877 energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
880 if (ir->efep != efepNO)
882 update_mdatoms(mdatoms, state->lambda[efptMASS]);
888 /* We need the kinetic energy at minus the half step for determining
889 * the full step kinetic energy and possibly for T-coupling.*/
890 /* This may not be quite working correctly yet . . . . */
891 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
892 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
893 nrnb, &vcm, wcycle, enerd, nullptr, nullptr, nullptr, nullptr, constr,
894 &nullSignaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
895 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
896 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
897 &top, makeConstArrayRef(state->x), state->box,
898 &shouldCheckNumberOfBondedInteractions);
900 clear_mat(force_vir);
902 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
904 /* Determine the energy and pressure:
905 * at nstcalcenergy steps and at energy output steps (set below).
907 if (EI_VV(ir->eI) && (!bInitStep))
909 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
910 bCalcVir = bCalcEnerStep
912 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
916 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
917 bCalcVir = bCalcEnerStep || (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
919 bCalcEner = bCalcEnerStep;
921 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
923 if (do_ene || do_log || bDoReplEx)
929 /* Do we need global communication ? */
930 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
931 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
933 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
934 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
935 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
939 /* Now is the time to relax the shells */
940 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
941 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
942 state->natoms, state->x.arrayRefWithPadding(),
943 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
944 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, graph,
945 shellfc, fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
949 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
950 is updated (or the AWH update will be performed twice for one step when continuing).
951 It would be best to call this update function from do_md_trajectory_writing but that
952 would occur after do_force. One would have to divide the update_awh function into one
953 function applying the AWH force and one doing the AWH bias update. The update AWH
954 bias function could then be called after do_md_trajectory_writing (then containing
955 update_awh_history). The checkpointing will in the future probably moved to the start
956 of the md loop which will rid of this issue. */
957 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
959 awh->updateHistory(state_global->awhHistory.get());
962 /* The coordinates (x) are shifted (to get whole molecules)
964 * This is parallellized as well, and does communication too.
965 * Check comments in sim_util.c
967 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation, imdSession, pull_work, step,
968 nrnb, wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
969 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, graph,
970 fr, runScheduleWork, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
971 (bNS ? GMX_FORCE_NS : 0) | force_flags, ddBalanceRegionHandler);
974 // VV integrators do not need the following velocity half step
975 // if it is the first step after starting from a checkpoint.
976 // That is, the half step is needed on all other steps, and
977 // also the first step when starting from a .tpr file.
978 if (EI_VV(ir->eI) && (!bFirstStep || startingBehavior == StartingBehavior::NewSimulation))
979 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
981 rvec* vbuf = nullptr;
983 wallcycle_start(wcycle, ewcUPDATE);
984 if (ir->eI == eiVV && bInitStep)
986 /* if using velocity verlet with full time step Ekin,
987 * take the first half step only to compute the
988 * virial for the first step. From there,
989 * revert back to the initial coordinates
990 * so that the input is actually the initial step.
992 snew(vbuf, state->natoms);
993 copy_rvecn(state->v.rvec_array(), vbuf, 0,
994 state->natoms); /* should make this better for parallelizing? */
998 /* this is for NHC in the Ekin(t+dt/2) version of vv */
999 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1000 trotter_seq, ettTSEQ1);
1003 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1004 etrtVELOCITY1, cr, constr);
1006 wallcycle_stop(wcycle, ewcUPDATE);
1007 constrain_velocities(step, nullptr, state, shake_vir, constr, bCalcVir, do_log, do_ene);
1008 wallcycle_start(wcycle, ewcUPDATE);
1009 /* if VV, compute the pressure and constraints */
1010 /* For VV2, we strictly only need this if using pressure
1011 * control, but we really would like to have accurate pressures
1013 * Think about ways around this in the future?
1014 * For now, keep this choice in comments.
1016 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1017 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1019 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1020 if (bCalcEner && ir->eI == eiVVAK)
1022 bSumEkinhOld = TRUE;
1024 /* for vv, the first half of the integration actually corresponds to the previous step.
1025 So we need information from the last step in the first half of the integration */
1026 if (bGStat || do_per_step(step - 1, nstglobalcomm))
1028 wallcycle_stop(wcycle, ewcUPDATE);
1030 gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1031 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
1032 nrnb, &vcm, wcycle, enerd, force_vir, shake_vir, total_vir, pres, constr,
1033 &nullSignaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1034 (bGStat ? CGLO_GSTAT : 0) | (bCalcEner ? CGLO_ENERGY : 0)
1035 | (bTemp ? CGLO_TEMPERATURE : 0) | (bPres ? CGLO_PRESSURE : 0)
1036 | (bPres ? CGLO_CONSTRAINT : 0) | (bStopCM ? CGLO_STOPCM : 0)
1037 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1040 /* explanation of above:
1041 a) We compute Ekin at the full time step
1042 if 1) we are using the AveVel Ekin, and it's not the
1043 initial step, or 2) if we are using AveEkin, but need the full
1044 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1045 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1046 EkinAveVel because it's needed for the pressure */
1047 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1048 top_global, &top, makeConstArrayRef(state->x),
1049 state->box, &shouldCheckNumberOfBondedInteractions);
1052 process_and_stopcm_grp(fplog, &vcm, *mdatoms, makeArrayRef(state->x),
1053 makeArrayRef(state->v));
1054 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1056 wallcycle_start(wcycle, ewcUPDATE);
1058 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1063 m_add(force_vir, shake_vir,
1064 total_vir); /* we need the un-dispersion corrected total vir here */
1065 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1066 trotter_seq, ettTSEQ2);
1068 /* TODO This is only needed when we're about to write
1069 * a checkpoint, because we use it after the restart
1070 * (in a kludge?). But what should we be doing if
1071 * the startingBehavior is NewSimulation or bInitStep are true? */
1072 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1074 copy_mat(shake_vir, state->svir_prev);
1075 copy_mat(force_vir, state->fvir_prev);
1077 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1079 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1080 enerd->term[F_TEMP] =
1081 sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1082 enerd->term[F_EKIN] = trace(ekind->ekin);
1085 else if (bExchanged)
1087 wallcycle_stop(wcycle, ewcUPDATE);
1088 /* We need the kinetic energy at minus the half step for determining
1089 * the full step kinetic energy and possibly for T-coupling.*/
1090 /* This may not be quite working correctly yet . . . . */
1091 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1092 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW],
1093 mdatoms, nrnb, &vcm, wcycle, enerd, nullptr, nullptr, nullptr,
1094 nullptr, constr, &nullSignaller, state->box, nullptr,
1095 &bSumEkinhOld, CGLO_GSTAT | CGLO_TEMPERATURE);
1096 wallcycle_start(wcycle, ewcUPDATE);
1099 /* if it's the initial step, we performed this first step just to get the constraint virial */
1100 if (ir->eI == eiVV && bInitStep)
1102 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1105 wallcycle_stop(wcycle, ewcUPDATE);
1108 /* compute the conserved quantity */
1111 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1114 last_ekin = enerd->term[F_EKIN];
1116 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1118 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1120 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1121 if (ir->efep != efepNO)
1123 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1127 /* ######## END FIRST UPDATE STEP ############## */
1128 /* ######## If doing VV, we now have v(dt) ###### */
1131 /* perform extended ensemble sampling in lambda - we don't
1132 actually move to the new state before outputting
1133 statistics, but if performing simulated tempering, we
1134 do update the velocities and the tau_t. */
1136 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state,
1137 state->dfhist, step, state->v.rvec_array(), mdatoms);
1138 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1141 copy_df_history(state_global->dfhist, state->dfhist);
1145 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1146 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1147 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1148 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1149 || checkpointHandler->isCheckpointingStep()))
1151 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1152 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1154 // Copy velocities if needed for the output/checkpointing.
1155 // NOTE: Copy on the search steps is done at the beginning of the step.
1156 if (useGpuForUpdate && !bNS
1157 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1159 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1160 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1162 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1163 // and update is offloaded hence forces are kept on the GPU for update and have not been
1164 // already transferred in do_force().
1165 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1166 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1167 // prior to GPU update.
1168 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1169 // copy call in do_force(...).
1170 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1171 // on host after the D2H copy in do_force(...).
1172 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1173 && do_per_step(step, ir->nstfout))
1175 stateGpu->copyForcesFromGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1176 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1178 /* Now we have the energies and forces corresponding to the
1179 * coordinates at time t. We must output all of this before
1182 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state, state_global,
1183 observablesHistory, top_global, fr, outf, energyOutput, ekind, f,
1184 checkpointHandler->isCheckpointingStep(), bRerunMD, bLastStep,
1185 mdrunOptions.writeConfout, bSumEkinhOld);
1186 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1187 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1189 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1190 if (startingBehavior != StartingBehavior::NewSimulation
1191 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1193 copy_mat(state->svir_prev, shake_vir);
1194 copy_mat(state->fvir_prev, force_vir);
1197 stopHandler->setSignal();
1198 resetHandler->setSignal(walltime_accounting);
1200 if (bGStat || !PAR(cr))
1202 /* In parallel we only have to check for checkpointing in steps
1203 * where we do global communication,
1204 * otherwise the other nodes don't know.
1206 checkpointHandler->setSignal(walltime_accounting);
1209 /* ######### START SECOND UPDATE STEP ################# */
1211 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1212 controlled in preprocessing */
1214 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1216 gmx_bool bIfRandomize;
1217 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1218 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1219 if (constr && bIfRandomize)
1221 constrain_velocities(step, nullptr, state, tmp_vir, constr, bCalcVir, do_log, do_ene);
1224 /* Box is changed in update() when we do pressure coupling,
1225 * but we should still use the old box for energy corrections and when
1226 * writing it to the energy file, so it matches the trajectory files for
1227 * the same timestep above. Make a copy in a separate array.
1229 copy_mat(state->box, lastbox);
1233 wallcycle_start(wcycle, ewcUPDATE);
1234 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1237 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1238 /* We can only do Berendsen coupling after we have summed
1239 * the kinetic energy or virial. Since the happens
1240 * in global_state after update, we should only do it at
1241 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1246 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1247 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1252 /* velocity half-step update */
1253 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1254 etrtVELOCITY2, cr, constr);
1257 /* Above, initialize just copies ekinh into ekin,
1258 * it doesn't copy position (for VV),
1259 * and entire integrator for MD.
1262 if (ir->eI == eiVVAK)
1264 cbuf.resize(state->x.size());
1265 std::copy(state->x.begin(), state->x.end(), cbuf.begin());
1268 /* With leap-frog type integrators we compute the kinetic energy
1269 * at a whole time step as the average of the half-time step kinetic
1270 * energies of two subsequent steps. Therefore we need to compute the
1271 * half step kinetic energy also if we need energies at the next step.
1273 const bool needHalfStepKineticEnergy =
1274 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1276 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1277 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1278 const bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN
1279 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1281 if (useGpuForUpdate)
1283 if (bNS && (bFirstStep || DOMAINDECOMP(cr)))
1285 integrator->set(stateGpu->getCoordinates(), stateGpu->getVelocities(),
1286 stateGpu->getForces(), top.idef, *mdatoms, ekind->ngtc);
1288 // Copy data to the GPU after buffers might have being reinitialized
1289 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1290 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1293 // If the buffer ops were not offloaded this step, the forces are on the host and have to be copied
1294 if (!runScheduleWork->stepWork.useGpuFBufferOps)
1296 stateGpu->copyForcesToGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1299 const bool doTemperatureScaling =
1300 (ir->etc != etcNO && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1302 // This applies Leap-Frog, LINCS and SETTLE in succession
1303 integrator->integrate(stateGpu->getForcesReadyOnDeviceEvent(
1304 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1305 ir->delta_t, true, bCalcVir, shake_vir, doTemperatureScaling,
1306 ekind->tcstat, doParrinelloRahman, ir->nstpcouple * ir->delta_t, M);
1308 // Copy velocities D2H after update if:
1309 // - Globals are computed this step (includes the energy output steps).
1310 // - Temperature is needed for the next step.
1311 if (bGStat || needHalfStepKineticEnergy)
1313 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1314 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1319 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1320 etrtPOSITION, cr, constr);
1322 wallcycle_stop(wcycle, ewcUPDATE);
1324 constrain_coordinates(step, &dvdl_constr, state, shake_vir, &upd, constr, bCalcVir,
1327 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state, cr, nrnb, wcycle, &upd,
1328 constr, do_log, do_ene);
1329 finish_update(ir, mdatoms, state, graph, nrnb, wcycle, &upd, constr);
1332 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1334 updatePrevStepPullCom(pull_work, state);
1337 if (ir->eI == eiVVAK)
1339 /* erase F_EKIN and F_TEMP here? */
1340 /* just compute the kinetic energy at the half step to perform a trotter step */
1341 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1342 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW],
1343 mdatoms, nrnb, &vcm, wcycle, enerd, force_vir, shake_vir, total_vir,
1344 pres, constr, &nullSignaller, lastbox, nullptr, &bSumEkinhOld,
1345 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE);
1346 wallcycle_start(wcycle, ewcUPDATE);
1347 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1348 /* now we know the scaling, we can compute the positions again */
1349 std::copy(cbuf.begin(), cbuf.end(), state->x.begin());
1351 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1352 etrtPOSITION, cr, constr);
1353 wallcycle_stop(wcycle, ewcUPDATE);
1355 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1356 /* are the small terms in the shake_vir here due
1357 * to numerical errors, or are they important
1358 * physically? I'm thinking they are just errors, but not completely sure.
1359 * For now, will call without actually constraining, constr=NULL*/
1360 finish_update(ir, mdatoms, state, graph, nrnb, wcycle, &upd, nullptr);
1364 /* this factor or 2 correction is necessary
1365 because half of the constraint force is removed
1366 in the vv step, so we have to double it. See
1367 the Redmine issue #1255. It is not yet clear
1368 if the factor of 2 is exact, or just a very
1369 good approximation, and this will be
1370 investigated. The next step is to see if this
1371 can be done adding a dhdl contribution from the
1372 rattle step, but this is somewhat more
1373 complicated with the current code. Will be
1374 investigated, hopefully for 4.6.3. However,
1375 this current solution is much better than
1376 having it completely wrong.
1378 enerd->term[F_DVDL_CONSTR] += 2 * dvdl_constr;
1382 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1385 if (vsite != nullptr)
1387 wallcycle_start(wcycle, ewcVSITECONSTR);
1388 if (graph != nullptr)
1390 shift_self(graph, state->box, state->x.rvec_array());
1392 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, state->v.rvec_array(),
1393 top.idef.iparams, top.idef.il, fr->pbcType, fr->bMolPBC, cr, state->box);
1395 if (graph != nullptr)
1397 unshift_self(graph, state->box, state->x.rvec_array());
1399 wallcycle_stop(wcycle, ewcVSITECONSTR);
1402 /* ############## IF NOT VV, Calculate globals HERE ############ */
1403 /* With Leap-Frog we can skip compute_globals at
1404 * non-communication steps, but we need to calculate
1405 * the kinetic energy one step before communication.
1408 // Organize to do inter-simulation signalling on steps if
1409 // and when algorithms require it.
1410 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1412 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1414 // Copy coordinates when needed to stop the CM motion.
1415 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1417 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1418 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1420 // Since we're already communicating at this step, we
1421 // can propagate intra-simulation signals. Note that
1422 // check_nstglobalcomm has the responsibility for
1423 // choosing the value of nstglobalcomm that is one way
1424 // bGStat becomes true, so we can't get into a
1425 // situation where e.g. checkpointing can't be
1427 bool doIntraSimSignal = true;
1428 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1431 gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1432 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
1433 nrnb, &vcm, wcycle, enerd, force_vir, shake_vir, total_vir, pres, constr,
1434 &signaller, lastbox, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1435 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1436 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1437 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1438 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1439 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1441 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1442 top_global, &top, makeConstArrayRef(state->x),
1443 state->box, &shouldCheckNumberOfBondedInteractions);
1444 if (!EI_VV(ir->eI) && bStopCM)
1446 process_and_stopcm_grp(fplog, &vcm, *mdatoms, makeArrayRef(state->x),
1447 makeArrayRef(state->v));
1448 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1450 // TODO: The special case of removing CM motion should be dealt more gracefully
1451 if (useGpuForUpdate)
1453 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1454 // Here we block until the H2D copy completes because event sync with the
1455 // force kernels that use the coordinates on the next steps is not implemented
1456 // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1457 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1463 /* ############# END CALC EKIN AND PRESSURE ################# */
1465 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1466 the virial that should probably be addressed eventually. state->veta has better properies,
1467 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1468 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1470 if (ir->efep != efepNO && !EI_VV(ir->eI))
1472 /* Sum up the foreign energy and dhdl terms for md and sd.
1473 Currently done every step so that dhdl is correct in the .edr */
1474 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1477 update_pcouple_after_coordinates(fplog, step, ir, mdatoms, pres, force_vir, shake_vir,
1478 pressureCouplingMu, state, nrnb, &upd, !useGpuForUpdate);
1480 const bool doBerendsenPressureCoupling =
1481 (inputrec->epc == epcBERENDSEN && do_per_step(step, inputrec->nstpcouple));
1482 if (useGpuForUpdate && (doBerendsenPressureCoupling || doParrinelloRahman))
1484 integrator->scaleCoordinates(pressureCouplingMu);
1485 integrator->setPbc(PbcType::Xyz, state->box);
1488 /* ################# END UPDATE STEP 2 ################# */
1489 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1491 /* The coordinates (x) were unshifted in update */
1494 /* We will not sum ekinh_old,
1495 * so signal that we still have to do it.
1497 bSumEkinhOld = TRUE;
1502 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1504 /* use the directly determined last velocity, not actually the averaged half steps */
1505 if (bTrotter && ir->eI == eiVV)
1507 enerd->term[F_EKIN] = last_ekin;
1509 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1511 if (integratorHasConservedEnergyQuantity(ir))
1515 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1519 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1522 /* ######### END PREPARING EDR OUTPUT ########### */
1528 if (fplog && do_log && bDoExpanded)
1530 /* only needed if doing expanded ensemble */
1531 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals,
1532 ir->bSimTemp ? ir->simtempvals : nullptr,
1533 state_global->dfhist, state->fep_state, ir->nstlog, step);
1537 energyOutput.addDataAtEnergyStep(bDoDHDL, bCalcEnerStep, t, mdatoms->tmass, enerd, state,
1538 ir->fepvals, ir->expandedvals, lastbox, shake_vir,
1539 force_vir, total_vir, pres, ekind, mu_tot, constr);
1543 energyOutput.recordNonEnergyStep();
1546 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1547 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1549 if (doSimulatedAnnealing)
1551 energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
1553 if (do_log || do_ene || do_dr || do_or)
1555 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1556 do_log ? fplog : nullptr, step, t, fcd, awh.get());
1561 pull_print_output(pull_work, step, t);
1564 if (do_per_step(step, ir->nstlog))
1566 if (fflush(fplog) != 0)
1568 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1574 /* Have to do this part _after_ outputting the logfile and the edr file */
1575 /* Gets written into the state at the beginning of next loop*/
1576 state->fep_state = lamnew;
1578 /* Print the remaining wall clock time for the run */
1579 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1583 fprintf(stderr, "\n");
1585 print_time(stderr, walltime_accounting, step, ir, cr);
1588 /* Ion/water position swapping.
1589 * Not done in last step since trajectory writing happens before this call
1590 * in the MD loop and exchanges would be lost anyway. */
1591 bNeedRepartition = FALSE;
1592 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep && do_per_step(step, ir->swap->nstswap))
1595 do_swapcoords(cr, step, t, ir, swap, wcycle, as_rvec_array(state->x.data()),
1596 state->box, MASTER(cr) && mdrunOptions.verbose, bRerunMD);
1598 if (bNeedRepartition && DOMAINDECOMP(cr))
1600 dd_collect_state(cr->dd, state, state_global);
1604 /* Replica exchange */
1608 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1611 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1613 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1, state_global, *top_global, ir,
1614 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
1615 nrnb, wcycle, FALSE);
1616 shouldCheckNumberOfBondedInteractions = true;
1617 upd.setNumAtoms(state->natoms);
1623 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1624 /* With all integrators, except VV, we need to retain the pressure
1625 * at the current step for coupling at the next step.
1627 if ((state->flags & (1U << estPRES_PREV))
1628 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1630 /* Store the pressure in t_state for pressure coupling
1631 * at the next MD step.
1633 copy_mat(pres, state->pres_prev);
1636 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1638 if ((membed != nullptr) && (!bLastStep))
1640 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1643 cycles = wallcycle_stop(wcycle, ewcSTEP);
1644 if (DOMAINDECOMP(cr) && wcycle)
1646 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1649 /* increase the MD step number */
1653 resetHandler->resetCounters(step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb,
1654 fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1656 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1657 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1659 /* End of main MD loop */
1661 /* Closing TNG files can include compressing data. Therefore it is good to do that
1662 * before stopping the time measurements. */
1663 mdoutf_tng_close(outf);
1665 /* Stop measuring walltime */
1666 walltime_accounting_end_time(walltime_accounting);
1668 if (!thisRankHasDuty(cr, DUTY_PME))
1670 /* Tell the PME only node to finish */
1671 gmx_pme_send_finish(cr);
1676 if (ir->nstcalcenergy > 0)
1678 energyOutput.printAnnealingTemperatures(fplog, groups, &(ir->opts));
1679 energyOutput.printAverages(fplog, groups);
1686 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1689 done_shellfc(fplog, shellfc, step_rel);
1691 if (useReplicaExchange && MASTER(cr))
1693 print_replica_exchange_statistics(fplog, repl_ex);
1696 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1698 global_stat_destroy(gstat);