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/pbc.h"
127 #include "gromacs/pulling/output.h"
128 #include "gromacs/pulling/pull.h"
129 #include "gromacs/swap/swapcoords.h"
130 #include "gromacs/timing/wallcycle.h"
131 #include "gromacs/timing/walltime_accounting.h"
132 #include "gromacs/topology/atoms.h"
133 #include "gromacs/topology/idef.h"
134 #include "gromacs/topology/mtop_util.h"
135 #include "gromacs/topology/topology.h"
136 #include "gromacs/trajectory/trajectoryframe.h"
137 #include "gromacs/utility/basedefinitions.h"
138 #include "gromacs/utility/cstringutil.h"
139 #include "gromacs/utility/fatalerror.h"
140 #include "gromacs/utility/logger.h"
141 #include "gromacs/utility/real.h"
142 #include "gromacs/utility/smalloc.h"
144 #include "legacysimulator.h"
145 #include "replicaexchange.h"
149 # include "corewrap.h"
152 using gmx::SimulationSignaller;
154 void gmx::LegacySimulator::do_md()
156 // TODO Historically, the EM and MD "integrators" used different
157 // names for the t_inputrec *parameter, but these must have the
158 // same name, now that it's a member of a struct. We use this ir
159 // alias to avoid a large ripple of nearly useless changes.
160 // t_inputrec is being replaced by IMdpOptionsProvider, so this
161 // will go away eventually.
162 t_inputrec* ir = inputrec;
163 int64_t step, step_rel;
164 double t, t0 = ir->init_t, lam0[efptNR];
165 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
166 gmx_bool bNS, bNStList, bStopCM, bFirstStep, bInitStep, bLastStep = FALSE;
167 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
168 gmx_bool do_ene, do_log, do_verbose;
169 gmx_bool bMasterState;
170 unsigned int force_flags;
171 tensor force_vir = { { 0 } }, shake_vir = { { 0 } }, total_vir = { { 0 } }, tmp_vir = { { 0 } },
175 matrix pressureCouplingMu, M;
176 gmx_repl_ex_t repl_ex = nullptr;
177 PaddedHostVector<gmx::RVec> f{};
178 gmx_global_stat_t gstat;
179 gmx_shellfc_t* shellfc;
180 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
181 gmx_bool bTemp, bPres, bTrotter;
183 std::vector<RVec> cbuf;
189 real saved_conserved_quantity = 0;
192 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
194 /* PME load balancing data for GPU kernels */
195 gmx_bool bPMETune = FALSE;
196 gmx_bool bPMETunePrinting = FALSE;
198 bool bInteractiveMDstep = false;
200 /* Domain decomposition could incorrectly miss a bonded
201 interaction, but checking for that requires a global
202 communication stage, which does not otherwise happen in DD
203 code. So we do that alongside the first global energy reduction
204 after a new DD is made. These variables handle whether the
205 check happens, and the result it returns. */
206 bool shouldCheckNumberOfBondedInteractions = false;
207 int totalNumberOfBondedInteractions = -1;
209 SimulationSignals signals;
210 // Most global communnication stages don't propagate mdrun
211 // signals, and will use this object to achieve that.
212 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
214 if (!mdrunOptions.writeConfout)
216 // This is on by default, and the main known use case for
217 // turning it off is for convenience in benchmarking, which is
218 // something that should not show up in the general user
223 "The -noconfout functionality is deprecated, and may be removed in a "
227 /* md-vv uses averaged full step velocities for T-control
228 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
229 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
230 bTrotter = (EI_VV(ir->eI)
231 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
233 const bool bRerunMD = false;
235 int nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
236 bGStatEveryStep = (nstglobalcomm == 1);
238 SimulationGroups* groups = &top_global->groups;
240 std::unique_ptr<EssentialDynamics> ed = nullptr;
241 if (opt2bSet("-ei", nfile, fnm))
243 /* Initialize essential dynamics sampling */
244 ed = init_edsam(mdlog, opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm), top_global,
245 ir, cr, constr, state_global, observablesHistory, oenv, startingBehavior);
247 else if (observablesHistory->edsamHistory)
250 "The checkpoint is from a run with essential dynamics sampling, "
251 "but the current run did not specify the -ei option. "
252 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
255 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
256 Update upd(ir, deform);
257 const bool doSimulatedAnnealing = initSimulatedAnnealing(ir, &upd);
258 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
260 bool simulationsShareState = false;
261 int nstSignalComm = nstglobalcomm;
263 // TODO This implementation of ensemble orientation restraints is nasty because
264 // a user can't just do multi-sim with single-sim orientation restraints.
265 bool usingEnsembleRestraints =
266 (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
267 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
269 // Replica exchange, ensemble restraints and AWH need all
270 // simulations to remain synchronized, so they need
271 // checkpoints and stop conditions to act on the same step, so
272 // the propagation of such signals must take place between
273 // simulations, not just within simulations.
274 // TODO: Make algorithm initializers set these flags.
275 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
277 if (simulationsShareState)
279 // Inter-simulation signal communication does not need to happen
280 // often, so we use a minimum of 200 steps to reduce overhead.
281 const int c_minimumInterSimulationSignallingInterval = 200;
282 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
287 if (startingBehavior != StartingBehavior::RestartWithAppending)
289 pleaseCiteCouplingAlgorithms(fplog, *ir);
292 init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, ir,
293 top_global, oenv, wcycle, startingBehavior, simulationsShareState, ms);
294 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
295 mdoutf_get_fp_dhdl(outf), false, startingBehavior, mdModulesNotifier);
297 gstat = global_stat_init(ir);
299 /* Check for polarizable models and flexible constraints */
300 shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
301 ir->nstcalcenergy, DOMAINDECOMP(cr));
304 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
305 if ((io > 2000) && MASTER(cr))
307 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
311 // Local state only becomes valid now.
312 std::unique_ptr<t_state> stateInstance;
316 gmx_localtop_t top(top_global->ffparams);
318 auto mdatoms = mdAtoms->mdatoms();
320 std::unique_ptr<UpdateConstrainGpu> integrator;
322 if (DOMAINDECOMP(cr))
324 stateInstance = std::make_unique<t_state>();
325 state = stateInstance.get();
326 dd_init_local_state(cr->dd, state_global, state);
328 /* Distribute the charge groups over the nodes from the master node */
329 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
330 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
331 nrnb, nullptr, FALSE);
332 shouldCheckNumberOfBondedInteractions = true;
333 upd.setNumAtoms(state->natoms);
337 state_change_natoms(state_global, state_global->natoms);
338 f.resizeWithPadding(state_global->natoms);
339 /* Copy the pointer to the global state */
340 state = state_global;
342 /* Generate and initialize new topology */
343 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, mdAtoms, constr, vsite, shellfc);
345 upd.setNumAtoms(state->natoms);
348 const auto& simulationWork = runScheduleWork->simulationWork;
349 const bool useGpuForPme = simulationWork.useGpuPme;
350 const bool useGpuForNonbonded = simulationWork.useGpuNonbonded;
351 const bool useGpuForBufferOps = simulationWork.useGpuBufferOps;
352 const bool useGpuForUpdate = simulationWork.useGpuUpdate;
354 StatePropagatorDataGpu* stateGpu = fr->stateGpu;
358 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
359 || constr->numConstraintsTotal() == 0,
360 "Constraints in domain decomposition are only supported with update "
361 "groups if using GPU update.\n");
362 GMX_RELEASE_ASSERT(ir->eConstrAlg != econtSHAKE || constr == nullptr
363 || constr->numConstraintsTotal() == 0,
364 "SHAKE is not supported with GPU update.");
365 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
366 "Either PME or short-ranged non-bonded interaction tasks must run on "
367 "the GPU to use GPU update.\n");
368 GMX_RELEASE_ASSERT(ir->eI == eiMD,
369 "Only the md integrator is supported with the GPU update.\n");
371 ir->etc != etcNOSEHOOVER,
372 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
373 GMX_RELEASE_ASSERT(ir->epc == epcNO || ir->epc == epcPARRINELLORAHMAN || ir->epc == epcBERENDSEN,
374 "Only Parrinello-Rahman and Berendsen pressure coupling are supported "
375 "with the GPU update.\n");
376 GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
377 "Virtual sites are not supported with the GPU update.\n");
378 GMX_RELEASE_ASSERT(ed == nullptr,
379 "Essential dynamics is not supported with the GPU update.\n");
380 GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(ir->pull),
381 "Constraints pulling is not supported with the GPU update.\n");
382 GMX_RELEASE_ASSERT(fcd->orires.nr == 0,
383 "Orientation restraints are not supported with the GPU update.\n");
386 || (!haveFreeEnergyType(*ir, efptBONDED) && !haveFreeEnergyType(*ir, efptMASS)),
387 "Free energy perturbation of masses and constraints are not supported with the GPU "
390 if (constr != nullptr && constr->numConstraintsTotal() > 0)
394 .appendText("Updating coordinates and applying constraints on the GPU.");
398 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
401 GMX_RELEASE_ASSERT(fr->deviceContext != nullptr,
402 "GPU device context should be initialized to use GPU update.");
403 GMX_RELEASE_ASSERT(stateGpu->getUpdateStream() != nullptr,
404 "Update stream can not be nullptr when update is on a GPU.");
405 integrator = std::make_unique<UpdateConstrainGpu>(*ir, *top_global, *fr->deviceContext,
406 *stateGpu->getUpdateStream(),
407 stateGpu->xUpdatedOnDevice());
409 integrator->setPbc(PbcType::Xyz, state->box);
412 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
414 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
416 if ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
418 changePinningPolicy(&f, PinningPolicy::PinnedIfSupported);
422 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
425 // NOTE: The global state is no longer used at this point.
426 // But state_global is still used as temporary storage space for writing
427 // the global state to file and potentially for replica exchange.
428 // (Global topology should persist.)
430 update_mdatoms(mdatoms, state->lambda[efptMASS]);
434 /* Check nstexpanded here, because the grompp check was broken */
435 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
438 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
440 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
445 EnergyElement::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
448 preparePrevStepPullCom(ir, pull_work, mdatoms, state, state_global, cr,
449 startingBehavior != StartingBehavior::NewSimulation);
451 // TODO: Remove this by converting AWH into a ForceProvider
452 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms,
453 startingBehavior != StartingBehavior::NewSimulation,
454 shellfc != nullptr, opt2fn("-awh", nfile, fnm), pull_work);
456 if (useReplicaExchange && MASTER(cr))
458 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir, replExParams);
460 /* PME tuning is only supported in the Verlet scheme, with PME for
461 * Coulomb. It is not supported with only LJ PME. */
462 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
463 && ir->cutoff_scheme != ecutsGROUP);
465 pme_load_balancing_t* pme_loadbal = nullptr;
468 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata,
472 if (!ir->bContinuation)
474 if (state->flags & (1U << estV))
476 auto v = makeArrayRef(state->v);
477 /* Set the velocities of vsites, shells and frozen atoms to zero */
478 for (i = 0; i < mdatoms->homenr; i++)
480 if (mdatoms->ptype[i] == eptVSite || mdatoms->ptype[i] == eptShell)
484 else if (mdatoms->cFREEZE)
486 for (m = 0; m < DIM; m++)
488 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
499 /* Constrain the initial coordinates and velocities */
500 do_constrain_first(fplog, constr, ir, mdatoms, state->natoms, state->x.arrayRefWithPadding(),
501 state->v.arrayRefWithPadding(), state->box, state->lambda[efptBONDED]);
505 /* Construct the virtual sites for the initial configuration */
506 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, nullptr, top.idef.iparams,
507 top.idef.il, fr->pbcType, fr->bMolPBC, cr, state->box);
511 if (ir->efep != efepNO)
513 /* Set free energy calculation frequency as the greatest common
514 * denominator of nstdhdl and repl_ex_nst. */
515 nstfep = ir->fepvals->nstdhdl;
518 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
520 if (useReplicaExchange)
522 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
526 /* Be REALLY careful about what flags you set here. You CANNOT assume
527 * this is the first step, since we might be restarting from a checkpoint,
528 * and in that case we should not do any modifications to the state.
530 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
532 // When restarting from a checkpoint, it can be appropriate to
533 // initialize ekind from quantities in the checkpoint. Otherwise,
534 // compute_globals must initialize ekind before the simulation
535 // starts/restarts. However, only the master rank knows what was
536 // found in the checkpoint file, so we have to communicate in
537 // order to coordinate the restart.
539 // TODO Consider removing this communication if/when checkpoint
540 // reading directly follows .tpr reading, because all ranks can
541 // agree on hasReadEkinState at that time.
542 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
545 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr);
547 if (hasReadEkinState)
549 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
552 unsigned int cglo_flags =
553 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
554 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
556 bSumEkinhOld = FALSE;
558 t_vcm vcm(top_global->groups, *ir);
559 reportComRemovalInfo(fplog, vcm);
561 /* To minimize communication, compute_globals computes the COM velocity
562 * and the kinetic energy for the velocities without COM motion removed.
563 * Thus to get the kinetic energy without the COM contribution, we need
564 * to call compute_globals twice.
566 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
568 unsigned int cglo_flags_iteration = cglo_flags;
569 if (bStopCM && cgloIteration == 0)
571 cglo_flags_iteration |= CGLO_STOPCM;
572 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
574 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
575 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
576 nrnb, &vcm, nullptr, enerd, force_vir, shake_vir, total_vir, pres, constr,
577 &nullSignaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
579 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
581 if (cglo_flags_iteration & CGLO_STOPCM)
583 /* At initialization, do not pass x with acceleration-correction mode
584 * to avoid (incorrect) correction of the initial coordinates.
586 auto x = (vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION) ? ArrayRef<RVec>()
587 : makeArrayRef(state->x);
588 process_and_stopcm_grp(fplog, &vcm, *mdatoms, x, makeArrayRef(state->v));
589 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
592 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
593 makeConstArrayRef(state->x), state->box,
594 &shouldCheckNumberOfBondedInteractions);
595 if (ir->eI == eiVVAK)
597 /* a second call to get the half step temperature initialized as well */
598 /* we do the same call as above, but turn the pressure off -- internally to
599 compute_globals, this is recognized as a velocity verlet half-step
600 kinetic energy calculation. This minimized excess variables, but
601 perhaps loses some logic?*/
603 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
604 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
605 nrnb, &vcm, nullptr, enerd, force_vir, shake_vir, total_vir, pres, constr,
606 &nullSignaller, state->box, nullptr, &bSumEkinhOld, cglo_flags & ~CGLO_PRESSURE);
609 /* Calculate the initial half step temperature, and save the ekinh_old */
610 if (startingBehavior == StartingBehavior::NewSimulation)
612 for (i = 0; (i < ir->opts.ngtc); i++)
614 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
618 /* need to make an initiation call to get the Trotter variables set, as well as other constants
619 for non-trotter temperature control */
620 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
624 if (!ir->bContinuation)
626 if (constr && ir->eConstrAlg == econtLINCS)
628 fprintf(fplog, "RMS relative constraint deviation after constraining: %.2e\n",
631 if (EI_STATE_VELOCITY(ir->eI))
633 real temp = enerd->term[F_TEMP];
636 /* Result of Ekin averaged over velocities of -half
637 * and +half step, while we only have -half step here.
641 fprintf(fplog, "Initial temperature: %g K\n", temp);
646 fprintf(stderr, "starting mdrun '%s'\n", *(top_global->name));
649 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
653 sprintf(tbuf, "%s", "infinite");
655 if (ir->init_step > 0)
657 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
658 gmx_step_str(ir->init_step + ir->nsteps, sbuf), tbuf,
659 gmx_step_str(ir->init_step, sbuf2), ir->init_step * ir->delta_t);
663 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
665 fprintf(fplog, "\n");
668 walltime_accounting_start_time(walltime_accounting);
669 wallcycle_start(wcycle, ewcRUN);
670 print_start(fplog, cr, walltime_accounting, "mdrun");
673 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
674 int chkpt_ret = fcCheckPointParallel(cr->nodeid, NULL, 0);
677 gmx_fatal(3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0);
681 /***********************************************************
685 ************************************************************/
688 /* Skip the first Nose-Hoover integration when we get the state from tpx */
689 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
690 bSumEkinhOld = FALSE;
692 bNeedRepartition = FALSE;
694 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
695 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
696 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
697 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
699 auto checkpointHandler = std::make_unique<CheckpointHandler>(
700 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]), simulationsShareState,
701 ir->nstlist == 0, MASTER(cr), mdrunOptions.writeConfout,
702 mdrunOptions.checkpointOptions.period);
704 const bool resetCountersIsLocal = true;
705 auto resetHandler = std::make_unique<ResetHandler>(
706 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
707 !resetCountersIsLocal, ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
708 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
710 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
712 step = ir->init_step;
715 // TODO extract this to new multi-simulation module
716 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
718 if (!multisim_int_all_are_equal(ms, ir->nsteps))
720 GMX_LOG(mdlog.warning)
722 "Note: The number of steps is not consistent across multi "
724 "but we are proceeding anyway!");
726 if (!multisim_int_all_are_equal(ms, ir->init_step))
728 if (simulationsShareState)
733 "The initial step is not consistent across multi simulations which "
740 GMX_LOG(mdlog.warning)
742 "Note: The initial step is not consistent across multi "
744 "but we are proceeding anyway!");
749 /* and stop now if we should */
750 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
754 /* Determine if this is a neighbor search step */
755 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
757 if (bPMETune && bNStList)
759 // This has to be here because PME load balancing is called so early.
760 // TODO: Move to after all booleans are defined.
761 if (useGpuForUpdate && !bFirstStep)
763 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
764 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
766 /* PME grid + cut-off optimization with GPUs or PME nodes */
767 pme_loadbal_do(pme_loadbal, cr, (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
768 fplog, mdlog, *ir, fr, state->box, state->x, wcycle, step, step_rel,
769 &bPMETunePrinting, simulationWork.useGpuPmePpCommunication);
772 wallcycle_start(wcycle, ewcSTEP);
774 bLastStep = (step_rel == ir->nsteps);
775 t = t0 + step * ir->delta_t;
777 // TODO Refactor this, so that nstfep does not need a default value of zero
778 if (ir->efep != efepNO || ir->bSimTemp)
780 /* find and set the current lambdas */
781 setCurrentLambdasLocal(step, ir->fepvals, lam0, state->lambda, state->fep_state);
783 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
784 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
785 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
786 && (step > 0) && (startingBehavior == StartingBehavior::NewSimulation));
789 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
790 && do_per_step(step, replExParams.exchangeInterval));
792 if (doSimulatedAnnealing)
794 update_annealing_target_temp(ir, t, &upd);
797 /* Stop Center of Mass motion */
798 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
800 /* Determine whether or not to do Neighbour Searching */
801 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
803 /* Note that the stopHandler will cause termination at nstglobalcomm
804 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
805 * nstpcouple steps, we have computed the half-step kinetic energy
806 * of the previous step and can always output energies at the last step.
808 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
810 /* do_log triggers energy and virial calculation. Because this leads
811 * to different code paths, forces can be different. Thus for exact
812 * continuation we should avoid extra log output.
813 * Note that the || bLastStep can result in non-exact continuation
814 * beyond the last step. But we don't consider that to be an issue.
816 do_log = (do_per_step(step, ir->nstlog)
817 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
818 do_verbose = mdrunOptions.verbose
819 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
821 if (useGpuForUpdate && !bFirstStep && bNS)
823 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
824 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
825 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
826 // Copy coordinate from the GPU when needed at the search step.
827 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
828 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
829 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
830 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
833 if (bNS && !(bFirstStep && ir->bContinuation))
835 bMasterState = FALSE;
836 /* Correct the new box if it is too skewed */
837 if (inputrecDynamicBox(ir))
839 if (correct_box(fplog, step, state->box))
842 // If update is offloaded, it should be informed about the box size change
845 integrator->setPbc(PbcType::Xyz, state->box);
849 if (DOMAINDECOMP(cr) && bMasterState)
851 dd_collect_state(cr->dd, state, state_global);
854 if (DOMAINDECOMP(cr))
856 /* Repartition the domain decomposition */
857 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
858 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
859 fr, vsite, constr, nrnb, wcycle, do_verbose && !bPMETunePrinting);
860 shouldCheckNumberOfBondedInteractions = true;
861 upd.setNumAtoms(state->natoms);
863 // Allocate or re-size GPU halo exchange object, if necessary
864 if (havePPDomainDecomposition(cr) && simulationWork.useGpuHaloExchange
865 && useGpuForNonbonded && is1D(*cr->dd))
867 // TODO remove need to pass local stream into GPU halo exchange - Redmine #3093
868 const DeviceStream* localStream =
869 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local);
870 const DeviceStream* nonLocalStream = Nbnxm::gpu_get_command_stream(
871 fr->nbv->gpu_nbv, InteractionLocality::NonLocal);
873 fr->deviceContext != nullptr,
874 "GPU device context should be initialized to use GPU halo exchange.");
875 GMX_RELEASE_ASSERT(localStream != nullptr,
876 "Local non-bonded stream can't be nullptr when using GPU "
878 GMX_RELEASE_ASSERT(nonLocalStream != nullptr,
879 "Non-local non-bonded stream can't be nullptr when using "
880 "GPU halo exchange.");
881 constructGpuHaloExchange(mdlog, *cr, *fr->deviceContext, *localStream, *nonLocalStream);
886 if (MASTER(cr) && do_log)
888 energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
891 if (ir->efep != efepNO)
893 update_mdatoms(mdatoms, state->lambda[efptMASS]);
899 /* We need the kinetic energy at minus the half step for determining
900 * the full step kinetic energy and possibly for T-coupling.*/
901 /* This may not be quite working correctly yet . . . . */
902 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
903 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
904 nrnb, &vcm, wcycle, enerd, nullptr, nullptr, nullptr, nullptr, constr,
905 &nullSignaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
906 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
907 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
908 &top, makeConstArrayRef(state->x), state->box,
909 &shouldCheckNumberOfBondedInteractions);
911 clear_mat(force_vir);
913 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
915 /* Determine the energy and pressure:
916 * at nstcalcenergy steps and at energy output steps (set below).
918 if (EI_VV(ir->eI) && (!bInitStep))
920 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
921 bCalcVir = bCalcEnerStep
923 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
927 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
928 bCalcVir = bCalcEnerStep || (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
930 bCalcEner = bCalcEnerStep;
932 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
934 if (do_ene || do_log || bDoReplEx)
940 /* Do we need global communication ? */
941 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
942 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
944 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
945 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
946 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
950 /* Now is the time to relax the shells */
951 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
952 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
953 state->natoms, state->x.arrayRefWithPadding(),
954 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
955 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, shellfc,
956 fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
960 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
961 is updated (or the AWH update will be performed twice for one step when continuing).
962 It would be best to call this update function from do_md_trajectory_writing but that
963 would occur after do_force. One would have to divide the update_awh function into one
964 function applying the AWH force and one doing the AWH bias update. The update AWH
965 bias function could then be called after do_md_trajectory_writing (then containing
966 update_awh_history). The checkpointing will in the future probably moved to the start
967 of the md loop which will rid of this issue. */
968 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
970 awh->updateHistory(state_global->awhHistory.get());
973 /* The coordinates (x) are shifted (to get whole molecules)
975 * This is parallellized as well, and does communication too.
976 * Check comments in sim_util.c
978 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation, imdSession, pull_work, step,
979 nrnb, wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
980 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, fr,
981 runScheduleWork, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
982 (bNS ? GMX_FORCE_NS : 0) | force_flags, ddBalanceRegionHandler);
985 // VV integrators do not need the following velocity half step
986 // if it is the first step after starting from a checkpoint.
987 // That is, the half step is needed on all other steps, and
988 // also the first step when starting from a .tpr file.
989 if (EI_VV(ir->eI) && (!bFirstStep || startingBehavior == StartingBehavior::NewSimulation))
990 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
992 rvec* vbuf = nullptr;
994 wallcycle_start(wcycle, ewcUPDATE);
995 if (ir->eI == eiVV && bInitStep)
997 /* if using velocity verlet with full time step Ekin,
998 * take the first half step only to compute the
999 * virial for the first step. From there,
1000 * revert back to the initial coordinates
1001 * so that the input is actually the initial step.
1003 snew(vbuf, state->natoms);
1004 copy_rvecn(state->v.rvec_array(), vbuf, 0,
1005 state->natoms); /* should make this better for parallelizing? */
1009 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1010 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1011 trotter_seq, ettTSEQ1);
1014 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1015 etrtVELOCITY1, cr, constr);
1017 wallcycle_stop(wcycle, ewcUPDATE);
1018 constrain_velocities(step, nullptr, state, shake_vir, constr, bCalcVir, do_log, do_ene);
1019 wallcycle_start(wcycle, ewcUPDATE);
1020 /* if VV, compute the pressure and constraints */
1021 /* For VV2, we strictly only need this if using pressure
1022 * control, but we really would like to have accurate pressures
1024 * Think about ways around this in the future?
1025 * For now, keep this choice in comments.
1027 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1028 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1030 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1031 if (bCalcEner && ir->eI == eiVVAK)
1033 bSumEkinhOld = TRUE;
1035 /* for vv, the first half of the integration actually corresponds to the previous step.
1036 So we need information from the last step in the first half of the integration */
1037 if (bGStat || do_per_step(step - 1, nstglobalcomm))
1039 wallcycle_stop(wcycle, ewcUPDATE);
1041 gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1042 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
1043 nrnb, &vcm, wcycle, enerd, force_vir, shake_vir, total_vir, pres, constr,
1044 &nullSignaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1045 (bGStat ? CGLO_GSTAT : 0) | (bCalcEner ? CGLO_ENERGY : 0)
1046 | (bTemp ? CGLO_TEMPERATURE : 0) | (bPres ? CGLO_PRESSURE : 0)
1047 | (bPres ? CGLO_CONSTRAINT : 0) | (bStopCM ? CGLO_STOPCM : 0)
1048 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1051 /* explanation of above:
1052 a) We compute Ekin at the full time step
1053 if 1) we are using the AveVel Ekin, and it's not the
1054 initial step, or 2) if we are using AveEkin, but need the full
1055 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1056 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1057 EkinAveVel because it's needed for the pressure */
1058 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1059 top_global, &top, makeConstArrayRef(state->x),
1060 state->box, &shouldCheckNumberOfBondedInteractions);
1063 process_and_stopcm_grp(fplog, &vcm, *mdatoms, makeArrayRef(state->x),
1064 makeArrayRef(state->v));
1065 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1067 wallcycle_start(wcycle, ewcUPDATE);
1069 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1074 m_add(force_vir, shake_vir,
1075 total_vir); /* we need the un-dispersion corrected total vir here */
1076 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1077 trotter_seq, ettTSEQ2);
1079 /* TODO This is only needed when we're about to write
1080 * a checkpoint, because we use it after the restart
1081 * (in a kludge?). But what should we be doing if
1082 * the startingBehavior is NewSimulation or bInitStep are true? */
1083 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1085 copy_mat(shake_vir, state->svir_prev);
1086 copy_mat(force_vir, state->fvir_prev);
1088 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1090 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1091 enerd->term[F_TEMP] =
1092 sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1093 enerd->term[F_EKIN] = trace(ekind->ekin);
1096 else if (bExchanged)
1098 wallcycle_stop(wcycle, ewcUPDATE);
1099 /* We need the kinetic energy at minus the half step for determining
1100 * the full step kinetic energy and possibly for T-coupling.*/
1101 /* This may not be quite working correctly yet . . . . */
1102 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1103 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW],
1104 mdatoms, nrnb, &vcm, wcycle, enerd, nullptr, nullptr, nullptr,
1105 nullptr, constr, &nullSignaller, state->box, nullptr,
1106 &bSumEkinhOld, CGLO_GSTAT | CGLO_TEMPERATURE);
1107 wallcycle_start(wcycle, ewcUPDATE);
1110 /* if it's the initial step, we performed this first step just to get the constraint virial */
1111 if (ir->eI == eiVV && bInitStep)
1113 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1116 wallcycle_stop(wcycle, ewcUPDATE);
1119 /* compute the conserved quantity */
1122 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1125 last_ekin = enerd->term[F_EKIN];
1127 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1129 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1131 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1132 if (ir->efep != efepNO)
1134 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1138 /* ######## END FIRST UPDATE STEP ############## */
1139 /* ######## If doing VV, we now have v(dt) ###### */
1142 /* perform extended ensemble sampling in lambda - we don't
1143 actually move to the new state before outputting
1144 statistics, but if performing simulated tempering, we
1145 do update the velocities and the tau_t. */
1147 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state,
1148 state->dfhist, step, state->v.rvec_array(), mdatoms);
1149 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1152 copy_df_history(state_global->dfhist, state->dfhist);
1156 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1157 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1158 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1159 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1160 || checkpointHandler->isCheckpointingStep()))
1162 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1163 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1165 // Copy velocities if needed for the output/checkpointing.
1166 // NOTE: Copy on the search steps is done at the beginning of the step.
1167 if (useGpuForUpdate && !bNS
1168 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1170 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1171 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1173 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1174 // and update is offloaded hence forces are kept on the GPU for update and have not been
1175 // already transferred in do_force().
1176 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1177 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1178 // prior to GPU update.
1179 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1180 // copy call in do_force(...).
1181 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1182 // on host after the D2H copy in do_force(...).
1183 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1184 && do_per_step(step, ir->nstfout))
1186 stateGpu->copyForcesFromGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1187 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1189 /* Now we have the energies and forces corresponding to the
1190 * coordinates at time t. We must output all of this before
1193 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state, state_global,
1194 observablesHistory, top_global, fr, outf, energyOutput, ekind, f,
1195 checkpointHandler->isCheckpointingStep(), bRerunMD, bLastStep,
1196 mdrunOptions.writeConfout, bSumEkinhOld);
1197 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1198 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1200 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1201 if (startingBehavior != StartingBehavior::NewSimulation
1202 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1204 copy_mat(state->svir_prev, shake_vir);
1205 copy_mat(state->fvir_prev, force_vir);
1208 stopHandler->setSignal();
1209 resetHandler->setSignal(walltime_accounting);
1211 if (bGStat || !PAR(cr))
1213 /* In parallel we only have to check for checkpointing in steps
1214 * where we do global communication,
1215 * otherwise the other nodes don't know.
1217 checkpointHandler->setSignal(walltime_accounting);
1220 /* ######### START SECOND UPDATE STEP ################# */
1222 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1223 controlled in preprocessing */
1225 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1227 gmx_bool bIfRandomize;
1228 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1229 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1230 if (constr && bIfRandomize)
1232 constrain_velocities(step, nullptr, state, tmp_vir, constr, bCalcVir, do_log, do_ene);
1235 /* Box is changed in update() when we do pressure coupling,
1236 * but we should still use the old box for energy corrections and when
1237 * writing it to the energy file, so it matches the trajectory files for
1238 * the same timestep above. Make a copy in a separate array.
1240 copy_mat(state->box, lastbox);
1244 wallcycle_start(wcycle, ewcUPDATE);
1245 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1248 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1249 /* We can only do Berendsen coupling after we have summed
1250 * the kinetic energy or virial. Since the happens
1251 * in global_state after update, we should only do it at
1252 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1257 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1258 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1263 /* velocity half-step update */
1264 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1265 etrtVELOCITY2, cr, constr);
1268 /* Above, initialize just copies ekinh into ekin,
1269 * it doesn't copy position (for VV),
1270 * and entire integrator for MD.
1273 if (ir->eI == eiVVAK)
1275 cbuf.resize(state->x.size());
1276 std::copy(state->x.begin(), state->x.end(), cbuf.begin());
1279 /* With leap-frog type integrators we compute the kinetic energy
1280 * at a whole time step as the average of the half-time step kinetic
1281 * energies of two subsequent steps. Therefore we need to compute the
1282 * half step kinetic energy also if we need energies at the next step.
1284 const bool needHalfStepKineticEnergy =
1285 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1287 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1288 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1289 const bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN
1290 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1292 if (useGpuForUpdate)
1294 if (bNS && (bFirstStep || DOMAINDECOMP(cr)))
1296 integrator->set(stateGpu->getCoordinates(), stateGpu->getVelocities(),
1297 stateGpu->getForces(), top.idef, *mdatoms, ekind->ngtc);
1299 // Copy data to the GPU after buffers might have being reinitialized
1300 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1301 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1304 // If the buffer ops were not offloaded this step, the forces are on the host and have to be copied
1305 if (!runScheduleWork->stepWork.useGpuFBufferOps)
1307 stateGpu->copyForcesToGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1310 const bool doTemperatureScaling =
1311 (ir->etc != etcNO && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1313 // This applies Leap-Frog, LINCS and SETTLE in succession
1314 integrator->integrate(stateGpu->getForcesReadyOnDeviceEvent(
1315 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1316 ir->delta_t, true, bCalcVir, shake_vir, doTemperatureScaling,
1317 ekind->tcstat, doParrinelloRahman, ir->nstpcouple * ir->delta_t, M);
1319 // Copy velocities D2H after update if:
1320 // - Globals are computed this step (includes the energy output steps).
1321 // - Temperature is needed for the next step.
1322 if (bGStat || needHalfStepKineticEnergy)
1324 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1325 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1330 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1331 etrtPOSITION, cr, constr);
1333 wallcycle_stop(wcycle, ewcUPDATE);
1335 constrain_coordinates(step, &dvdl_constr, state, shake_vir, &upd, constr, bCalcVir,
1338 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state, cr, nrnb, wcycle, &upd,
1339 constr, do_log, do_ene);
1340 finish_update(ir, mdatoms, state, wcycle, &upd, constr);
1343 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1345 updatePrevStepPullCom(pull_work, state);
1348 if (ir->eI == eiVVAK)
1350 /* erase F_EKIN and F_TEMP here? */
1351 /* just compute the kinetic energy at the half step to perform a trotter step */
1352 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1353 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW],
1354 mdatoms, nrnb, &vcm, wcycle, enerd, force_vir, shake_vir, total_vir,
1355 pres, constr, &nullSignaller, lastbox, nullptr, &bSumEkinhOld,
1356 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE);
1357 wallcycle_start(wcycle, ewcUPDATE);
1358 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1359 /* now we know the scaling, we can compute the positions again */
1360 std::copy(cbuf.begin(), cbuf.end(), state->x.begin());
1362 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1363 etrtPOSITION, cr, constr);
1364 wallcycle_stop(wcycle, ewcUPDATE);
1366 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1367 /* are the small terms in the shake_vir here due
1368 * to numerical errors, or are they important
1369 * physically? I'm thinking they are just errors, but not completely sure.
1370 * For now, will call without actually constraining, constr=NULL*/
1371 finish_update(ir, mdatoms, state, wcycle, &upd, nullptr);
1375 /* this factor or 2 correction is necessary
1376 because half of the constraint force is removed
1377 in the vv step, so we have to double it. See
1378 the Redmine issue #1255. It is not yet clear
1379 if the factor of 2 is exact, or just a very
1380 good approximation, and this will be
1381 investigated. The next step is to see if this
1382 can be done adding a dhdl contribution from the
1383 rattle step, but this is somewhat more
1384 complicated with the current code. Will be
1385 investigated, hopefully for 4.6.3. However,
1386 this current solution is much better than
1387 having it completely wrong.
1389 enerd->term[F_DVDL_CONSTR] += 2 * dvdl_constr;
1393 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1396 if (vsite != nullptr)
1398 wallcycle_start(wcycle, ewcVSITECONSTR);
1399 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, state->v.rvec_array(),
1400 top.idef.iparams, top.idef.il, fr->pbcType, fr->bMolPBC, cr, state->box);
1401 wallcycle_stop(wcycle, ewcVSITECONSTR);
1404 /* ############## IF NOT VV, Calculate globals HERE ############ */
1405 /* With Leap-Frog we can skip compute_globals at
1406 * non-communication steps, but we need to calculate
1407 * the kinetic energy one step before communication.
1410 // Organize to do inter-simulation signalling on steps if
1411 // and when algorithms require it.
1412 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1414 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1416 // Copy coordinates when needed to stop the CM motion.
1417 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1419 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1420 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1422 // Since we're already communicating at this step, we
1423 // can propagate intra-simulation signals. Note that
1424 // check_nstglobalcomm has the responsibility for
1425 // choosing the value of nstglobalcomm that is one way
1426 // bGStat becomes true, so we can't get into a
1427 // situation where e.g. checkpointing can't be
1429 bool doIntraSimSignal = true;
1430 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1433 gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1434 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
1435 nrnb, &vcm, wcycle, enerd, force_vir, shake_vir, total_vir, pres, constr,
1436 &signaller, lastbox, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1437 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1438 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1439 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1440 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1441 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1443 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1444 top_global, &top, makeConstArrayRef(state->x),
1445 state->box, &shouldCheckNumberOfBondedInteractions);
1446 if (!EI_VV(ir->eI) && bStopCM)
1448 process_and_stopcm_grp(fplog, &vcm, *mdatoms, makeArrayRef(state->x),
1449 makeArrayRef(state->v));
1450 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1452 // TODO: The special case of removing CM motion should be dealt more gracefully
1453 if (useGpuForUpdate)
1455 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1456 // Here we block until the H2D copy completes because event sync with the
1457 // force kernels that use the coordinates on the next steps is not implemented
1458 // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1459 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1465 /* ############# END CALC EKIN AND PRESSURE ################# */
1467 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1468 the virial that should probably be addressed eventually. state->veta has better properies,
1469 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1470 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1472 if (ir->efep != efepNO && !EI_VV(ir->eI))
1474 /* Sum up the foreign energy and dhdl terms for md and sd.
1475 Currently done every step so that dhdl is correct in the .edr */
1476 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1479 update_pcouple_after_coordinates(fplog, step, ir, mdatoms, pres, force_vir, shake_vir,
1480 pressureCouplingMu, state, nrnb, &upd, !useGpuForUpdate);
1482 const bool doBerendsenPressureCoupling =
1483 (inputrec->epc == epcBERENDSEN && do_per_step(step, inputrec->nstpcouple));
1484 if (useGpuForUpdate && (doBerendsenPressureCoupling || doParrinelloRahman))
1486 integrator->scaleCoordinates(pressureCouplingMu);
1487 integrator->setPbc(PbcType::Xyz, state->box);
1490 /* ################# END UPDATE STEP 2 ################# */
1491 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1493 /* The coordinates (x) were unshifted in update */
1496 /* We will not sum ekinh_old,
1497 * so signal that we still have to do it.
1499 bSumEkinhOld = TRUE;
1504 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1506 /* use the directly determined last velocity, not actually the averaged half steps */
1507 if (bTrotter && ir->eI == eiVV)
1509 enerd->term[F_EKIN] = last_ekin;
1511 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1513 if (integratorHasConservedEnergyQuantity(ir))
1517 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1521 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1524 /* ######### END PREPARING EDR OUTPUT ########### */
1530 if (fplog && do_log && bDoExpanded)
1532 /* only needed if doing expanded ensemble */
1533 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals,
1534 ir->bSimTemp ? ir->simtempvals : nullptr,
1535 state_global->dfhist, state->fep_state, ir->nstlog, step);
1539 energyOutput.addDataAtEnergyStep(bDoDHDL, bCalcEnerStep, t, mdatoms->tmass, enerd, state,
1540 ir->fepvals, ir->expandedvals, lastbox, shake_vir,
1541 force_vir, total_vir, pres, ekind, mu_tot, constr);
1545 energyOutput.recordNonEnergyStep();
1548 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1549 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1551 if (doSimulatedAnnealing)
1553 energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
1555 if (do_log || do_ene || do_dr || do_or)
1557 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1558 do_log ? fplog : nullptr, step, t, fcd, awh.get());
1563 pull_print_output(pull_work, step, t);
1566 if (do_per_step(step, ir->nstlog))
1568 if (fflush(fplog) != 0)
1570 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1576 /* Have to do this part _after_ outputting the logfile and the edr file */
1577 /* Gets written into the state at the beginning of next loop*/
1578 state->fep_state = lamnew;
1580 /* Print the remaining wall clock time for the run */
1581 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1585 fprintf(stderr, "\n");
1587 print_time(stderr, walltime_accounting, step, ir, cr);
1590 /* Ion/water position swapping.
1591 * Not done in last step since trajectory writing happens before this call
1592 * in the MD loop and exchanges would be lost anyway. */
1593 bNeedRepartition = FALSE;
1594 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep && do_per_step(step, ir->swap->nstswap))
1597 do_swapcoords(cr, step, t, ir, swap, wcycle, as_rvec_array(state->x.data()),
1598 state->box, MASTER(cr) && mdrunOptions.verbose, bRerunMD);
1600 if (bNeedRepartition && DOMAINDECOMP(cr))
1602 dd_collect_state(cr->dd, state, state_global);
1606 /* Replica exchange */
1610 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1613 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1615 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1, state_global, *top_global, ir,
1616 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
1617 nrnb, wcycle, FALSE);
1618 shouldCheckNumberOfBondedInteractions = true;
1619 upd.setNumAtoms(state->natoms);
1625 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1626 /* With all integrators, except VV, we need to retain the pressure
1627 * at the current step for coupling at the next step.
1629 if ((state->flags & (1U << estPRES_PREV))
1630 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1632 /* Store the pressure in t_state for pressure coupling
1633 * at the next MD step.
1635 copy_mat(pres, state->pres_prev);
1638 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1640 if ((membed != nullptr) && (!bLastStep))
1642 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1645 cycles = wallcycle_stop(wcycle, ewcSTEP);
1646 if (DOMAINDECOMP(cr) && wcycle)
1648 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1651 /* increase the MD step number */
1655 resetHandler->resetCounters(step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb,
1656 fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1658 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1659 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1661 /* End of main MD loop */
1663 /* Closing TNG files can include compressing data. Therefore it is good to do that
1664 * before stopping the time measurements. */
1665 mdoutf_tng_close(outf);
1667 /* Stop measuring walltime */
1668 walltime_accounting_end_time(walltime_accounting);
1670 if (!thisRankHasDuty(cr, DUTY_PME))
1672 /* Tell the PME only node to finish */
1673 gmx_pme_send_finish(cr);
1678 if (ir->nstcalcenergy > 0)
1680 energyOutput.printAnnealingTemperatures(fplog, groups, &(ir->opts));
1681 energyOutput.printAverages(fplog, groups);
1688 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1691 done_shellfc(fplog, shellfc, step_rel);
1693 if (useReplicaExchange && MASTER(cr))
1695 print_replica_exchange_statistics(fplog, repl_ex);
1698 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1700 global_stat_destroy(gstat);