Merge branch 'origin/release-2020' into master
[gromacs.git] / src / gromacs / mdrun / md.cpp
blob10b76df3778955f5a83038e18ce2edd367911be7
1 /*
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.
37 /*! \internal \file
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
44 #include "gmxpre.h"
46 #include <cinttypes>
47 #include <cmath>
48 #include <cstdio>
49 #include <cstdlib>
51 #include <algorithm>
52 #include <memory>
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/device_stream_manager.h"
71 #include "gromacs/gpu_utils/gpu_utils.h"
72 #include "gromacs/imd/imd.h"
73 #include "gromacs/listed_forces/manage_threading.h"
74 #include "gromacs/math/functions.h"
75 #include "gromacs/math/utilities.h"
76 #include "gromacs/math/vec.h"
77 #include "gromacs/math/vectypes.h"
78 #include "gromacs/mdlib/checkpointhandler.h"
79 #include "gromacs/mdlib/compute_io.h"
80 #include "gromacs/mdlib/constr.h"
81 #include "gromacs/mdlib/ebin.h"
82 #include "gromacs/mdlib/enerdata_utils.h"
83 #include "gromacs/mdlib/energyoutput.h"
84 #include "gromacs/mdlib/expanded.h"
85 #include "gromacs/mdlib/force.h"
86 #include "gromacs/mdlib/force_flags.h"
87 #include "gromacs/mdlib/forcerec.h"
88 #include "gromacs/mdlib/md_support.h"
89 #include "gromacs/mdlib/mdatoms.h"
90 #include "gromacs/mdlib/mdoutf.h"
91 #include "gromacs/mdlib/membed.h"
92 #include "gromacs/mdlib/resethandler.h"
93 #include "gromacs/mdlib/sighandler.h"
94 #include "gromacs/mdlib/simulationsignal.h"
95 #include "gromacs/mdlib/stat.h"
96 #include "gromacs/mdlib/stophandler.h"
97 #include "gromacs/mdlib/tgroup.h"
98 #include "gromacs/mdlib/trajectory_writing.h"
99 #include "gromacs/mdlib/update.h"
100 #include "gromacs/mdlib/update_constrain_gpu.h"
101 #include "gromacs/mdlib/vcm.h"
102 #include "gromacs/mdlib/vsite.h"
103 #include "gromacs/mdrunutility/handlerestart.h"
104 #include "gromacs/mdrunutility/multisim.h"
105 #include "gromacs/mdrunutility/printtime.h"
106 #include "gromacs/mdtypes/awh_history.h"
107 #include "gromacs/mdtypes/awh_params.h"
108 #include "gromacs/mdtypes/commrec.h"
109 #include "gromacs/mdtypes/df_history.h"
110 #include "gromacs/mdtypes/energyhistory.h"
111 #include "gromacs/mdtypes/fcdata.h"
112 #include "gromacs/mdtypes/forcerec.h"
113 #include "gromacs/mdtypes/group.h"
114 #include "gromacs/mdtypes/inputrec.h"
115 #include "gromacs/mdtypes/interaction_const.h"
116 #include "gromacs/mdtypes/md_enums.h"
117 #include "gromacs/mdtypes/mdatom.h"
118 #include "gromacs/mdtypes/mdrunoptions.h"
119 #include "gromacs/mdtypes/observableshistory.h"
120 #include "gromacs/mdtypes/pullhistory.h"
121 #include "gromacs/mdtypes/simulation_workload.h"
122 #include "gromacs/mdtypes/state.h"
123 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
124 #include "gromacs/modularsimulator/energyelement.h"
125 #include "gromacs/nbnxm/gpu_data_mgmt.h"
126 #include "gromacs/nbnxm/nbnxm.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"
147 #include "shellfc.h"
149 #if GMX_FAHCORE
150 # include "corewrap.h"
151 #endif
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 } },
173 pres = { { 0 } };
174 int i, m;
175 rvec mu_tot;
176 matrix pressureCouplingMu, M;
177 gmx_repl_ex_t repl_ex = nullptr;
178 PaddedHostVector<gmx::RVec> f{};
179 gmx_global_stat_t gstat;
180 gmx_shellfc_t* shellfc;
181 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
182 gmx_bool bTemp, bPres, bTrotter;
183 real dvdl_constr;
184 std::vector<RVec> cbuf;
185 matrix lastbox;
186 int lamnew = 0;
187 /* for FEP */
188 int nstfep = 0;
189 double cycles;
190 real saved_conserved_quantity = 0;
191 real last_ekin = 0;
192 t_extmass MassQ;
193 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
195 /* PME load balancing data for GPU kernels */
196 gmx_bool bPMETune = FALSE;
197 gmx_bool bPMETunePrinting = FALSE;
199 bool bInteractiveMDstep = false;
201 /* Domain decomposition could incorrectly miss a bonded
202 interaction, but checking for that requires a global
203 communication stage, which does not otherwise happen in DD
204 code. So we do that alongside the first global energy reduction
205 after a new DD is made. These variables handle whether the
206 check happens, and the result it returns. */
207 bool shouldCheckNumberOfBondedInteractions = false;
208 int totalNumberOfBondedInteractions = -1;
210 SimulationSignals signals;
211 // Most global communnication stages don't propagate mdrun
212 // signals, and will use this object to achieve that.
213 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
215 if (!mdrunOptions.writeConfout)
217 // This is on by default, and the main known use case for
218 // turning it off is for convenience in benchmarking, which is
219 // something that should not show up in the general user
220 // interface.
221 GMX_LOG(mdlog.info)
222 .asParagraph()
223 .appendText(
224 "The -noconfout functionality is deprecated, and may be removed in a "
225 "future version.");
228 /* md-vv uses averaged full step velocities for T-control
229 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
230 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
231 bTrotter = (EI_VV(ir->eI)
232 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
234 const bool bRerunMD = false;
236 int nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
237 bGStatEveryStep = (nstglobalcomm == 1);
239 const SimulationGroups* groups = &top_global->groups;
241 std::unique_ptr<EssentialDynamics> ed = nullptr;
242 if (opt2bSet("-ei", nfile, fnm))
244 /* Initialize essential dynamics sampling */
245 ed = init_edsam(mdlog, opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm), top_global,
246 ir, cr, constr, state_global, observablesHistory, oenv, startingBehavior);
248 else if (observablesHistory->edsamHistory)
250 gmx_fatal(FARGS,
251 "The checkpoint is from a run with essential dynamics sampling, "
252 "but the current run did not specify the -ei option. "
253 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
256 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
257 Update upd(ir, deform);
258 const bool doSimulatedAnnealing = initSimulatedAnnealing(ir, &upd);
259 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
261 bool simulationsShareState = false;
262 int nstSignalComm = nstglobalcomm;
264 // TODO This implementation of ensemble orientation restraints is nasty because
265 // a user can't just do multi-sim with single-sim orientation restraints.
266 bool usingEnsembleRestraints =
267 (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
268 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
270 // Replica exchange, ensemble restraints and AWH need all
271 // simulations to remain synchronized, so they need
272 // checkpoints and stop conditions to act on the same step, so
273 // the propagation of such signals must take place between
274 // simulations, not just within simulations.
275 // TODO: Make algorithm initializers set these flags.
276 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
278 if (simulationsShareState)
280 // Inter-simulation signal communication does not need to happen
281 // often, so we use a minimum of 200 steps to reduce overhead.
282 const int c_minimumInterSimulationSignallingInterval = 200;
283 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
284 * nstglobalcomm;
288 if (startingBehavior != StartingBehavior::RestartWithAppending)
290 pleaseCiteCouplingAlgorithms(fplog, *ir);
292 gmx_mdoutf* outf =
293 init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, ir,
294 top_global, oenv, wcycle, startingBehavior, simulationsShareState, ms);
295 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
296 mdoutf_get_fp_dhdl(outf), false, startingBehavior, mdModulesNotifier);
298 gstat = global_stat_init(ir);
300 /* Check for polarizable models and flexible constraints */
301 shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
302 ir->nstcalcenergy, DOMAINDECOMP(cr));
305 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
306 if ((io > 2000) && MASTER(cr))
308 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
312 // Local state only becomes valid now.
313 std::unique_ptr<t_state> stateInstance;
314 t_state* state;
317 gmx_localtop_t top(top_global->ffparams);
319 auto mdatoms = mdAtoms->mdatoms();
321 std::unique_ptr<UpdateConstrainGpu> integrator;
323 if (DOMAINDECOMP(cr))
325 stateInstance = std::make_unique<t_state>();
326 state = stateInstance.get();
327 dd_init_local_state(cr->dd, state_global, state);
329 /* Distribute the charge groups over the nodes from the master node */
330 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
331 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
332 nrnb, nullptr, FALSE);
333 shouldCheckNumberOfBondedInteractions = true;
334 upd.setNumAtoms(state->natoms);
336 else
338 state_change_natoms(state_global, 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, &f, 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;
356 // TODO: the assertions below should be handled by UpdateConstraintsBuilder.
357 if (useGpuForUpdate)
359 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
360 || constr->numConstraintsTotal() == 0,
361 "Constraints in domain decomposition are only supported with update "
362 "groups if using GPU update.\n");
363 GMX_RELEASE_ASSERT(ir->eConstrAlg != econtSHAKE || constr == nullptr
364 || constr->numConstraintsTotal() == 0,
365 "SHAKE is not supported with GPU update.");
366 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
367 "Either PME or short-ranged non-bonded interaction tasks must run on "
368 "the GPU to use GPU update.\n");
369 GMX_RELEASE_ASSERT(ir->eI == eiMD,
370 "Only the md integrator is supported with the GPU update.\n");
371 GMX_RELEASE_ASSERT(
372 ir->etc != etcNOSEHOOVER,
373 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
374 GMX_RELEASE_ASSERT(ir->epc == epcNO || ir->epc == epcPARRINELLORAHMAN || ir->epc == epcBERENDSEN,
375 "Only Parrinello-Rahman and Berendsen pressure coupling are supported "
376 "with the GPU update.\n");
377 GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
378 "Virtual sites are not supported with the GPU update.\n");
379 GMX_RELEASE_ASSERT(ed == nullptr,
380 "Essential dynamics is not supported with the GPU update.\n");
381 GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(ir->pull),
382 "Constraints pulling is not supported with the GPU update.\n");
383 GMX_RELEASE_ASSERT(fcd->orires.nr == 0,
384 "Orientation restraints are not supported with the GPU update.\n");
385 GMX_RELEASE_ASSERT(
386 ir->efep == efepNO
387 || (!haveFreeEnergyType(*ir, efptBONDED) && !haveFreeEnergyType(*ir, efptMASS)),
388 "Free energy perturbation of masses and constraints are not supported with the GPU "
389 "update.");
391 if (constr != nullptr && constr->numConstraintsTotal() > 0)
393 GMX_LOG(mdlog.info)
394 .asParagraph()
395 .appendText("Updating coordinates and applying constraints on the GPU.");
397 else
399 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
401 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
402 "Device stream manager should be initialized in order to use GPU "
403 "update-constraints.");
404 GMX_RELEASE_ASSERT(
405 fr->deviceStreamManager->streamIsValid(gmx::DeviceStreamType::UpdateAndConstraints),
406 "Update stream should be initialized in order to use GPU "
407 "update-constraints.");
408 integrator = std::make_unique<UpdateConstrainGpu>(
409 *ir, *top_global, fr->deviceStreamManager->context(),
410 fr->deviceStreamManager->stream(gmx::DeviceStreamType::UpdateAndConstraints),
411 stateGpu->xUpdatedOnDevice());
413 integrator->setPbc(PbcType::Xyz, state->box);
416 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
418 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
420 if ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
422 changePinningPolicy(&f, PinningPolicy::PinnedIfSupported);
424 if (useGpuForUpdate)
426 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
429 // NOTE: The global state is no longer used at this point.
430 // But state_global is still used as temporary storage space for writing
431 // the global state to file and potentially for replica exchange.
432 // (Global topology should persist.)
434 update_mdatoms(mdatoms, state->lambda[efptMASS]);
436 if (ir->bExpanded)
438 /* Check nstexpanded here, because the grompp check was broken */
439 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
441 gmx_fatal(FARGS,
442 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
444 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
447 if (MASTER(cr))
449 EnergyElement::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
452 preparePrevStepPullCom(ir, pull_work, mdatoms, state, state_global, cr,
453 startingBehavior != StartingBehavior::NewSimulation);
455 // TODO: Remove this by converting AWH into a ForceProvider
456 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms,
457 startingBehavior != StartingBehavior::NewSimulation,
458 shellfc != nullptr, opt2fn("-awh", nfile, fnm), pull_work);
460 if (useReplicaExchange && MASTER(cr))
462 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir, replExParams);
464 /* PME tuning is only supported in the Verlet scheme, with PME for
465 * Coulomb. It is not supported with only LJ PME. */
466 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
467 && ir->cutoff_scheme != ecutsGROUP);
469 pme_load_balancing_t* pme_loadbal = nullptr;
470 if (bPMETune)
472 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata,
473 fr->nbv->useGpu());
476 if (!ir->bContinuation)
478 if (state->flags & (1U << estV))
480 auto v = makeArrayRef(state->v);
481 /* Set the velocities of vsites, shells and frozen atoms to zero */
482 for (i = 0; i < mdatoms->homenr; i++)
484 if (mdatoms->ptype[i] == eptVSite || mdatoms->ptype[i] == eptShell)
486 clear_rvec(v[i]);
488 else if (mdatoms->cFREEZE)
490 for (m = 0; m < DIM; m++)
492 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
494 v[i][m] = 0;
501 if (constr)
503 /* Constrain the initial coordinates and velocities */
504 do_constrain_first(fplog, constr, ir, mdatoms, state->natoms, state->x.arrayRefWithPadding(),
505 state->v.arrayRefWithPadding(), state->box, state->lambda[efptBONDED]);
507 if (vsite)
509 /* Construct the virtual sites for the initial configuration */
510 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, nullptr, top.idef.iparams,
511 top.idef.il, fr->pbcType, fr->bMolPBC, cr, state->box);
515 if (ir->efep != efepNO)
517 /* Set free energy calculation frequency as the greatest common
518 * denominator of nstdhdl and repl_ex_nst. */
519 nstfep = ir->fepvals->nstdhdl;
520 if (ir->bExpanded)
522 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
524 if (useReplicaExchange)
526 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
530 /* Be REALLY careful about what flags you set here. You CANNOT assume
531 * this is the first step, since we might be restarting from a checkpoint,
532 * and in that case we should not do any modifications to the state.
534 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
536 // When restarting from a checkpoint, it can be appropriate to
537 // initialize ekind from quantities in the checkpoint. Otherwise,
538 // compute_globals must initialize ekind before the simulation
539 // starts/restarts. However, only the master rank knows what was
540 // found in the checkpoint file, so we have to communicate in
541 // order to coordinate the restart.
543 // TODO Consider removing this communication if/when checkpoint
544 // reading directly follows .tpr reading, because all ranks can
545 // agree on hasReadEkinState at that time.
546 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
547 if (PAR(cr))
549 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr->mpi_comm_mygroup);
551 if (hasReadEkinState)
553 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
556 unsigned int cglo_flags =
557 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
558 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
560 bSumEkinhOld = FALSE;
562 t_vcm vcm(top_global->groups, *ir);
563 reportComRemovalInfo(fplog, vcm);
565 /* To minimize communication, compute_globals computes the COM velocity
566 * and the kinetic energy for the velocities without COM motion removed.
567 * Thus to get the kinetic energy without the COM contribution, we need
568 * to call compute_globals twice.
570 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
572 unsigned int cglo_flags_iteration = cglo_flags;
573 if (bStopCM && cgloIteration == 0)
575 cglo_flags_iteration |= CGLO_STOPCM;
576 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
578 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
579 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
580 nrnb, &vcm, nullptr, enerd, force_vir, shake_vir, total_vir, pres, constr,
581 &nullSignaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
582 cglo_flags_iteration
583 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
584 : 0));
585 if (cglo_flags_iteration & CGLO_STOPCM)
587 /* At initialization, do not pass x with acceleration-correction mode
588 * to avoid (incorrect) correction of the initial coordinates.
590 auto x = (vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION) ? ArrayRef<RVec>()
591 : makeArrayRef(state->x);
592 process_and_stopcm_grp(fplog, &vcm, *mdatoms, x, makeArrayRef(state->v));
593 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
596 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
597 makeConstArrayRef(state->x), state->box,
598 &shouldCheckNumberOfBondedInteractions);
599 if (ir->eI == eiVVAK)
601 /* a second call to get the half step temperature initialized as well */
602 /* we do the same call as above, but turn the pressure off -- internally to
603 compute_globals, this is recognized as a velocity verlet half-step
604 kinetic energy calculation. This minimized excess variables, but
605 perhaps loses some logic?*/
607 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
608 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
609 nrnb, &vcm, nullptr, enerd, force_vir, shake_vir, total_vir, pres, constr,
610 &nullSignaller, state->box, nullptr, &bSumEkinhOld, cglo_flags & ~CGLO_PRESSURE);
613 /* Calculate the initial half step temperature, and save the ekinh_old */
614 if (startingBehavior == StartingBehavior::NewSimulation)
616 for (i = 0; (i < ir->opts.ngtc); i++)
618 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
622 /* need to make an initiation call to get the Trotter variables set, as well as other constants
623 for non-trotter temperature control */
624 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
626 if (MASTER(cr))
628 if (!ir->bContinuation)
630 if (constr && ir->eConstrAlg == econtLINCS)
632 fprintf(fplog, "RMS relative constraint deviation after constraining: %.2e\n",
633 constr->rmsd());
635 if (EI_STATE_VELOCITY(ir->eI))
637 real temp = enerd->term[F_TEMP];
638 if (ir->eI != eiVV)
640 /* Result of Ekin averaged over velocities of -half
641 * and +half step, while we only have -half step here.
643 temp *= 2;
645 fprintf(fplog, "Initial temperature: %g K\n", temp);
649 char tbuf[20];
650 fprintf(stderr, "starting mdrun '%s'\n", *(top_global->name));
651 if (ir->nsteps >= 0)
653 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
655 else
657 sprintf(tbuf, "%s", "infinite");
659 if (ir->init_step > 0)
661 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
662 gmx_step_str(ir->init_step + ir->nsteps, sbuf), tbuf,
663 gmx_step_str(ir->init_step, sbuf2), ir->init_step * ir->delta_t);
665 else
667 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
669 fprintf(fplog, "\n");
672 walltime_accounting_start_time(walltime_accounting);
673 wallcycle_start(wcycle, ewcRUN);
674 print_start(fplog, cr, walltime_accounting, "mdrun");
676 #if GMX_FAHCORE
677 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
678 int chkpt_ret = fcCheckPointParallel(cr->nodeid, NULL, 0);
679 if (chkpt_ret == 0)
681 gmx_fatal(3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0);
683 #endif
685 /***********************************************************
687 * Loop over MD steps
689 ************************************************************/
691 bFirstStep = TRUE;
692 /* Skip the first Nose-Hoover integration when we get the state from tpx */
693 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
694 bSumEkinhOld = FALSE;
695 bExchanged = FALSE;
696 bNeedRepartition = FALSE;
698 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
699 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
700 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
701 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
703 auto checkpointHandler = std::make_unique<CheckpointHandler>(
704 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]), simulationsShareState,
705 ir->nstlist == 0, MASTER(cr), mdrunOptions.writeConfout,
706 mdrunOptions.checkpointOptions.period);
708 const bool resetCountersIsLocal = true;
709 auto resetHandler = std::make_unique<ResetHandler>(
710 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
711 !resetCountersIsLocal, ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
712 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
714 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
716 step = ir->init_step;
717 step_rel = 0;
719 // TODO extract this to new multi-simulation module
720 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
722 if (!multisim_int_all_are_equal(ms, ir->nsteps))
724 GMX_LOG(mdlog.warning)
725 .appendText(
726 "Note: The number of steps is not consistent across multi "
727 "simulations,\n"
728 "but we are proceeding anyway!");
730 if (!multisim_int_all_are_equal(ms, ir->init_step))
732 if (simulationsShareState)
734 if (MASTER(cr))
736 gmx_fatal(FARGS,
737 "The initial step is not consistent across multi simulations which "
738 "share the state");
740 gmx_barrier(cr->mpi_comm_mygroup);
742 else
744 GMX_LOG(mdlog.warning)
745 .appendText(
746 "Note: The initial step is not consistent across multi "
747 "simulations,\n"
748 "but we are proceeding anyway!");
753 /* and stop now if we should */
754 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
755 while (!bLastStep)
758 /* Determine if this is a neighbor search step */
759 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
761 if (bPMETune && bNStList)
763 // This has to be here because PME load balancing is called so early.
764 // TODO: Move to after all booleans are defined.
765 if (useGpuForUpdate && !bFirstStep)
767 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
768 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
770 /* PME grid + cut-off optimization with GPUs or PME nodes */
771 pme_loadbal_do(pme_loadbal, cr, (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
772 fplog, mdlog, *ir, fr, state->box, state->x, wcycle, step, step_rel,
773 &bPMETunePrinting, simulationWork.useGpuPmePpCommunication);
776 wallcycle_start(wcycle, ewcSTEP);
778 bLastStep = (step_rel == ir->nsteps);
779 t = t0 + step * ir->delta_t;
781 // TODO Refactor this, so that nstfep does not need a default value of zero
782 if (ir->efep != efepNO || ir->bSimTemp)
784 /* find and set the current lambdas */
785 setCurrentLambdasLocal(step, ir->fepvals, lam0, state->lambda, state->fep_state);
787 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
788 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
789 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
790 && (!bFirstStep));
793 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
794 && do_per_step(step, replExParams.exchangeInterval));
796 if (doSimulatedAnnealing)
798 update_annealing_target_temp(ir, t, &upd);
801 /* Stop Center of Mass motion */
802 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
804 /* Determine whether or not to do Neighbour Searching */
805 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
807 /* Note that the stopHandler will cause termination at nstglobalcomm
808 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
809 * nstpcouple steps, we have computed the half-step kinetic energy
810 * of the previous step and can always output energies at the last step.
812 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
814 /* do_log triggers energy and virial calculation. Because this leads
815 * to different code paths, forces can be different. Thus for exact
816 * continuation we should avoid extra log output.
817 * Note that the || bLastStep can result in non-exact continuation
818 * beyond the last step. But we don't consider that to be an issue.
820 do_log = (do_per_step(step, ir->nstlog)
821 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
822 do_verbose = mdrunOptions.verbose
823 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
825 if (useGpuForUpdate && !bFirstStep && bNS)
827 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
828 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
829 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
830 // Copy coordinate from the GPU when needed at the search step.
831 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
832 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
833 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
834 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
837 if (bNS && !(bFirstStep && ir->bContinuation))
839 bMasterState = FALSE;
840 /* Correct the new box if it is too skewed */
841 if (inputrecDynamicBox(ir))
843 if (correct_box(fplog, step, state->box))
845 bMasterState = TRUE;
846 // If update is offloaded, it should be informed about the box size change
847 if (useGpuForUpdate)
849 integrator->setPbc(PbcType::Xyz, state->box);
853 if (DOMAINDECOMP(cr) && bMasterState)
855 dd_collect_state(cr->dd, state, state_global);
858 if (DOMAINDECOMP(cr))
860 /* Repartition the domain decomposition */
861 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
862 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
863 fr, vsite, constr, nrnb, wcycle, do_verbose && !bPMETunePrinting);
864 shouldCheckNumberOfBondedInteractions = true;
865 upd.setNumAtoms(state->natoms);
867 // Allocate or re-size GPU halo exchange object, if necessary
868 if (havePPDomainDecomposition(cr) && simulationWork.useGpuHaloExchange
869 && useGpuForNonbonded && is1D(*cr->dd))
871 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
872 "GPU device manager has to be initialized to use GPU "
873 "version of halo exchange.");
874 // TODO remove need to pass local stream into GPU halo exchange - Issue #3093
875 constructGpuHaloExchange(mdlog, *cr, *fr->deviceStreamManager);
880 if (MASTER(cr) && do_log)
882 gmx::EnergyOutput::printHeader(fplog, step,
883 t); /* can we improve the information printed here? */
886 if (ir->efep != efepNO)
888 update_mdatoms(mdatoms, state->lambda[efptMASS]);
891 if (bExchanged)
894 /* We need the kinetic energy at minus the half step for determining
895 * the full step kinetic energy and possibly for T-coupling.*/
896 /* This may not be quite working correctly yet . . . . */
897 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
898 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
899 nrnb, &vcm, wcycle, enerd, nullptr, nullptr, nullptr, nullptr, constr,
900 &nullSignaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
901 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
902 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
903 &top, makeConstArrayRef(state->x), state->box,
904 &shouldCheckNumberOfBondedInteractions);
906 clear_mat(force_vir);
908 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
910 /* Determine the energy and pressure:
911 * at nstcalcenergy steps and at energy output steps (set below).
913 if (EI_VV(ir->eI) && (!bInitStep))
915 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
916 bCalcVir = bCalcEnerStep
917 || (ir->epc != epcNO
918 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
920 else
922 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
923 bCalcVir = bCalcEnerStep || (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
925 bCalcEner = bCalcEnerStep;
927 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
929 if (do_ene || do_log || bDoReplEx)
931 bCalcVir = TRUE;
932 bCalcEner = TRUE;
935 /* Do we need global communication ? */
936 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
937 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
939 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
940 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
941 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
943 if (shellfc)
945 /* Now is the time to relax the shells */
946 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
947 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
948 state->natoms, state->x.arrayRefWithPadding(),
949 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
950 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, shellfc,
951 fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
953 else
955 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
956 is updated (or the AWH update will be performed twice for one step when continuing).
957 It would be best to call this update function from do_md_trajectory_writing but that
958 would occur after do_force. One would have to divide the update_awh function into one
959 function applying the AWH force and one doing the AWH bias update. The update AWH
960 bias function could then be called after do_md_trajectory_writing (then containing
961 update_awh_history). The checkpointing will in the future probably moved to the start
962 of the md loop which will rid of this issue. */
963 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
965 awh->updateHistory(state_global->awhHistory.get());
968 /* The coordinates (x) are shifted (to get whole molecules)
969 * in do_force.
970 * This is parallellized as well, and does communication too.
971 * Check comments in sim_util.c
973 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation, imdSession, pull_work, step,
974 nrnb, wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
975 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, fr,
976 runScheduleWork, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
977 (bNS ? GMX_FORCE_NS : 0) | force_flags, ddBalanceRegionHandler);
980 // VV integrators do not need the following velocity half step
981 // if it is the first step after starting from a checkpoint.
982 // That is, the half step is needed on all other steps, and
983 // also the first step when starting from a .tpr file.
984 if (EI_VV(ir->eI) && (!bFirstStep || startingBehavior == StartingBehavior::NewSimulation))
985 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
987 rvec* vbuf = nullptr;
989 wallcycle_start(wcycle, ewcUPDATE);
990 if (ir->eI == eiVV && bInitStep)
992 /* if using velocity verlet with full time step Ekin,
993 * take the first half step only to compute the
994 * virial for the first step. From there,
995 * revert back to the initial coordinates
996 * so that the input is actually the initial step.
998 snew(vbuf, state->natoms);
999 copy_rvecn(state->v.rvec_array(), vbuf, 0,
1000 state->natoms); /* should make this better for parallelizing? */
1002 else
1004 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1005 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1006 trotter_seq, ettTSEQ1);
1009 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1010 etrtVELOCITY1, cr, constr);
1012 wallcycle_stop(wcycle, ewcUPDATE);
1013 constrain_velocities(step, nullptr, state, shake_vir, constr, bCalcVir, do_log, do_ene);
1014 wallcycle_start(wcycle, ewcUPDATE);
1015 /* if VV, compute the pressure and constraints */
1016 /* For VV2, we strictly only need this if using pressure
1017 * control, but we really would like to have accurate pressures
1018 * printed out.
1019 * Think about ways around this in the future?
1020 * For now, keep this choice in comments.
1022 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1023 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1024 bPres = TRUE;
1025 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1026 if (bCalcEner && ir->eI == eiVVAK)
1028 bSumEkinhOld = TRUE;
1030 /* for vv, the first half of the integration actually corresponds to the previous step.
1031 So we need information from the last step in the first half of the integration */
1032 if (bGStat || do_per_step(step - 1, nstglobalcomm))
1034 wallcycle_stop(wcycle, ewcUPDATE);
1035 compute_globals(
1036 gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1037 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
1038 nrnb, &vcm, wcycle, enerd, force_vir, shake_vir, total_vir, pres, constr,
1039 &nullSignaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1040 (bGStat ? CGLO_GSTAT : 0) | (bCalcEner ? CGLO_ENERGY : 0)
1041 | (bTemp ? CGLO_TEMPERATURE : 0) | (bPres ? CGLO_PRESSURE : 0)
1042 | (bPres ? CGLO_CONSTRAINT : 0) | (bStopCM ? CGLO_STOPCM : 0)
1043 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1044 : 0)
1045 | CGLO_SCALEEKIN);
1046 /* explanation of above:
1047 a) We compute Ekin at the full time step
1048 if 1) we are using the AveVel Ekin, and it's not the
1049 initial step, or 2) if we are using AveEkin, but need the full
1050 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1051 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1052 EkinAveVel because it's needed for the pressure */
1053 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1054 top_global, &top, makeConstArrayRef(state->x),
1055 state->box, &shouldCheckNumberOfBondedInteractions);
1056 if (bStopCM)
1058 process_and_stopcm_grp(fplog, &vcm, *mdatoms, makeArrayRef(state->x),
1059 makeArrayRef(state->v));
1060 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1062 wallcycle_start(wcycle, ewcUPDATE);
1064 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1065 if (!bInitStep)
1067 if (bTrotter)
1069 m_add(force_vir, shake_vir,
1070 total_vir); /* we need the un-dispersion corrected total vir here */
1071 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1072 trotter_seq, ettTSEQ2);
1074 /* TODO This is only needed when we're about to write
1075 * a checkpoint, because we use it after the restart
1076 * (in a kludge?). But what should we be doing if
1077 * the startingBehavior is NewSimulation or bInitStep are true? */
1078 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1080 copy_mat(shake_vir, state->svir_prev);
1081 copy_mat(force_vir, state->fvir_prev);
1083 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1085 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1086 enerd->term[F_TEMP] =
1087 sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1088 enerd->term[F_EKIN] = trace(ekind->ekin);
1091 else if (bExchanged)
1093 wallcycle_stop(wcycle, ewcUPDATE);
1094 /* We need the kinetic energy at minus the half step for determining
1095 * the full step kinetic energy and possibly for T-coupling.*/
1096 /* This may not be quite working correctly yet . . . . */
1097 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1098 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW],
1099 mdatoms, nrnb, &vcm, wcycle, enerd, nullptr, nullptr, nullptr,
1100 nullptr, constr, &nullSignaller, state->box, nullptr,
1101 &bSumEkinhOld, CGLO_GSTAT | CGLO_TEMPERATURE);
1102 wallcycle_start(wcycle, ewcUPDATE);
1105 /* if it's the initial step, we performed this first step just to get the constraint virial */
1106 if (ir->eI == eiVV && bInitStep)
1108 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1109 sfree(vbuf);
1111 wallcycle_stop(wcycle, ewcUPDATE);
1114 /* compute the conserved quantity */
1115 if (EI_VV(ir->eI))
1117 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1118 if (ir->eI == eiVV)
1120 last_ekin = enerd->term[F_EKIN];
1122 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1124 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1126 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1127 if (ir->efep != efepNO)
1129 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1133 /* ######## END FIRST UPDATE STEP ############## */
1134 /* ######## If doing VV, we now have v(dt) ###### */
1135 if (bDoExpanded)
1137 /* perform extended ensemble sampling in lambda - we don't
1138 actually move to the new state before outputting
1139 statistics, but if performing simulated tempering, we
1140 do update the velocities and the tau_t. */
1142 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state,
1143 state->dfhist, step, state->v.rvec_array(), mdatoms);
1144 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1145 if (MASTER(cr))
1147 copy_df_history(state_global->dfhist, state->dfhist);
1151 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1152 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1153 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1154 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1155 || checkpointHandler->isCheckpointingStep()))
1157 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1158 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1160 // Copy velocities if needed for the output/checkpointing.
1161 // NOTE: Copy on the search steps is done at the beginning of the step.
1162 if (useGpuForUpdate && !bNS
1163 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1165 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1166 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1168 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1169 // and update is offloaded hence forces are kept on the GPU for update and have not been
1170 // already transferred in do_force().
1171 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1172 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1173 // prior to GPU update.
1174 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1175 // copy call in do_force(...).
1176 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1177 // on host after the D2H copy in do_force(...).
1178 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1179 && do_per_step(step, ir->nstfout))
1181 stateGpu->copyForcesFromGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1182 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1184 /* Now we have the energies and forces corresponding to the
1185 * coordinates at time t. We must output all of this before
1186 * the update.
1188 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state, state_global,
1189 observablesHistory, top_global, fr, outf, energyOutput, ekind, f,
1190 checkpointHandler->isCheckpointingStep(), bRerunMD, bLastStep,
1191 mdrunOptions.writeConfout, bSumEkinhOld);
1192 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1193 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1195 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1196 if (startingBehavior != StartingBehavior::NewSimulation && bFirstStep
1197 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1199 copy_mat(state->svir_prev, shake_vir);
1200 copy_mat(state->fvir_prev, force_vir);
1203 stopHandler->setSignal();
1204 resetHandler->setSignal(walltime_accounting);
1206 if (bGStat || !PAR(cr))
1208 /* In parallel we only have to check for checkpointing in steps
1209 * where we do global communication,
1210 * otherwise the other nodes don't know.
1212 checkpointHandler->setSignal(walltime_accounting);
1215 /* ######### START SECOND UPDATE STEP ################# */
1217 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1218 controlled in preprocessing */
1220 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1222 gmx_bool bIfRandomize;
1223 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1224 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1225 if (constr && bIfRandomize)
1227 constrain_velocities(step, nullptr, state, tmp_vir, constr, bCalcVir, do_log, do_ene);
1230 /* Box is changed in update() when we do pressure coupling,
1231 * but we should still use the old box for energy corrections and when
1232 * writing it to the energy file, so it matches the trajectory files for
1233 * the same timestep above. Make a copy in a separate array.
1235 copy_mat(state->box, lastbox);
1237 dvdl_constr = 0;
1239 wallcycle_start(wcycle, ewcUPDATE);
1240 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1241 if (bTrotter)
1243 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1244 /* We can only do Berendsen coupling after we have summed
1245 * the kinetic energy or virial. Since the happens
1246 * in global_state after update, we should only do it at
1247 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1250 else
1252 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1253 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1256 if (EI_VV(ir->eI))
1258 /* velocity half-step update */
1259 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1260 etrtVELOCITY2, cr, constr);
1263 /* Above, initialize just copies ekinh into ekin,
1264 * it doesn't copy position (for VV),
1265 * and entire integrator for MD.
1268 if (ir->eI == eiVVAK)
1270 cbuf.resize(state->x.size());
1271 std::copy(state->x.begin(), state->x.end(), cbuf.begin());
1274 /* With leap-frog type integrators we compute the kinetic energy
1275 * at a whole time step as the average of the half-time step kinetic
1276 * energies of two subsequent steps. Therefore we need to compute the
1277 * half step kinetic energy also if we need energies at the next step.
1279 const bool needHalfStepKineticEnergy =
1280 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1282 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1283 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1284 const bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN
1285 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1287 if (useGpuForUpdate)
1289 if (bNS && (bFirstStep || DOMAINDECOMP(cr)))
1291 integrator->set(stateGpu->getCoordinates(), stateGpu->getVelocities(),
1292 stateGpu->getForces(), top.idef, *mdatoms, ekind->ngtc);
1294 // Copy data to the GPU after buffers might have being reinitialized
1295 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1296 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1299 // If the buffer ops were not offloaded this step, the forces are on the host and have to be copied
1300 if (!runScheduleWork->stepWork.useGpuFBufferOps)
1302 stateGpu->copyForcesToGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1305 const bool doTemperatureScaling =
1306 (ir->etc != etcNO && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1308 // This applies Leap-Frog, LINCS and SETTLE in succession
1309 integrator->integrate(stateGpu->getForcesReadyOnDeviceEvent(
1310 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1311 ir->delta_t, true, bCalcVir, shake_vir, doTemperatureScaling,
1312 ekind->tcstat, doParrinelloRahman, ir->nstpcouple * ir->delta_t, M);
1314 // Copy velocities D2H after update if:
1315 // - Globals are computed this step (includes the energy output steps).
1316 // - Temperature is needed for the next step.
1317 if (bGStat || needHalfStepKineticEnergy)
1319 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1320 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1323 else
1325 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1326 etrtPOSITION, cr, constr);
1328 wallcycle_stop(wcycle, ewcUPDATE);
1330 constrain_coordinates(step, &dvdl_constr, state, shake_vir, &upd, constr, bCalcVir,
1331 do_log, do_ene);
1333 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state, cr, nrnb, wcycle, &upd,
1334 constr, do_log, do_ene);
1335 finish_update(ir, mdatoms, state, wcycle, &upd, constr);
1338 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1340 updatePrevStepPullCom(pull_work, state);
1343 if (ir->eI == eiVVAK)
1345 /* erase F_EKIN and F_TEMP here? */
1346 /* just compute the kinetic energy at the half step to perform a trotter step */
1347 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1348 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW],
1349 mdatoms, nrnb, &vcm, wcycle, enerd, force_vir, shake_vir, total_vir,
1350 pres, constr, &nullSignaller, lastbox, nullptr, &bSumEkinhOld,
1351 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE);
1352 wallcycle_start(wcycle, ewcUPDATE);
1353 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1354 /* now we know the scaling, we can compute the positions again */
1355 std::copy(cbuf.begin(), cbuf.end(), state->x.begin());
1357 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1358 etrtPOSITION, cr, constr);
1359 wallcycle_stop(wcycle, ewcUPDATE);
1361 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1362 /* are the small terms in the shake_vir here due
1363 * to numerical errors, or are they important
1364 * physically? I'm thinking they are just errors, but not completely sure.
1365 * For now, will call without actually constraining, constr=NULL*/
1366 finish_update(ir, mdatoms, state, wcycle, &upd, nullptr);
1368 if (EI_VV(ir->eI))
1370 /* this factor or 2 correction is necessary
1371 because half of the constraint force is removed
1372 in the vv step, so we have to double it. See
1373 the Issue #1255. It is not yet clear
1374 if the factor of 2 is exact, or just a very
1375 good approximation, and this will be
1376 investigated. The next step is to see if this
1377 can be done adding a dhdl contribution from the
1378 rattle step, but this is somewhat more
1379 complicated with the current code. Will be
1380 investigated, hopefully for 4.6.3. However,
1381 this current solution is much better than
1382 having it completely wrong.
1384 enerd->term[F_DVDL_CONSTR] += 2 * dvdl_constr;
1386 else
1388 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1391 if (vsite != nullptr)
1393 wallcycle_start(wcycle, ewcVSITECONSTR);
1394 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, state->v.rvec_array(),
1395 top.idef.iparams, top.idef.il, fr->pbcType, fr->bMolPBC, cr, state->box);
1396 wallcycle_stop(wcycle, ewcVSITECONSTR);
1399 /* ############## IF NOT VV, Calculate globals HERE ############ */
1400 /* With Leap-Frog we can skip compute_globals at
1401 * non-communication steps, but we need to calculate
1402 * the kinetic energy one step before communication.
1405 // Organize to do inter-simulation signalling on steps if
1406 // and when algorithms require it.
1407 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1409 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1411 // Copy coordinates when needed to stop the CM motion.
1412 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1414 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1415 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1417 // Since we're already communicating at this step, we
1418 // can propagate intra-simulation signals. Note that
1419 // check_nstglobalcomm has the responsibility for
1420 // choosing the value of nstglobalcomm that is one way
1421 // bGStat becomes true, so we can't get into a
1422 // situation where e.g. checkpointing can't be
1423 // signalled.
1424 bool doIntraSimSignal = true;
1425 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1427 compute_globals(
1428 gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1429 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
1430 nrnb, &vcm, wcycle, enerd, force_vir, shake_vir, total_vir, pres, constr,
1431 &signaller, lastbox, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1432 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1433 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1434 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1435 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1436 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1437 : 0));
1438 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1439 top_global, &top, makeConstArrayRef(state->x),
1440 state->box, &shouldCheckNumberOfBondedInteractions);
1441 if (!EI_VV(ir->eI) && bStopCM)
1443 process_and_stopcm_grp(fplog, &vcm, *mdatoms, makeArrayRef(state->x),
1444 makeArrayRef(state->v));
1445 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1447 // TODO: The special case of removing CM motion should be dealt more gracefully
1448 if (useGpuForUpdate)
1450 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1451 // Here we block until the H2D copy completes because event sync with the
1452 // force kernels that use the coordinates on the next steps is not implemented
1453 // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1454 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1455 // If the COM removal changed the velocities on the CPU, this has to be accounted for.
1456 if (vcm.mode != ecmNO)
1458 stateGpu->copyVelocitiesToGpu(state->v, 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 */
1494 if (!bGStat)
1496 /* We will not sum ekinh_old,
1497 * so signal that we still have to do it.
1499 bSumEkinhOld = TRUE;
1502 if (bCalcEner)
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))
1515 if (EI_VV(ir->eI))
1517 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1519 else
1521 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1524 /* ######### END PREPARING EDR OUTPUT ########### */
1527 /* Output stuff */
1528 if (MASTER(cr))
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);
1537 if (bCalcEner)
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);
1543 else
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 gmx::EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups,
1554 &(ir->opts));
1556 if (do_log || do_ene || do_dr || do_or)
1558 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1559 do_log ? fplog : nullptr, step, t, fcd, awh.get());
1562 if (ir->bPull)
1564 pull_print_output(pull_work, step, t);
1567 if (do_per_step(step, ir->nstlog))
1569 if (fflush(fplog) != 0)
1571 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1575 if (bDoExpanded)
1577 /* Have to do this part _after_ outputting the logfile and the edr file */
1578 /* Gets written into the state at the beginning of next loop*/
1579 state->fep_state = lamnew;
1581 /* Print the remaining wall clock time for the run */
1582 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1584 if (shellfc)
1586 fprintf(stderr, "\n");
1588 print_time(stderr, walltime_accounting, step, ir, cr);
1591 /* Ion/water position swapping.
1592 * Not done in last step since trajectory writing happens before this call
1593 * in the MD loop and exchanges would be lost anyway. */
1594 bNeedRepartition = FALSE;
1595 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep && do_per_step(step, ir->swap->nstswap))
1597 bNeedRepartition =
1598 do_swapcoords(cr, step, t, ir, swap, wcycle, as_rvec_array(state->x.data()),
1599 state->box, MASTER(cr) && mdrunOptions.verbose, bRerunMD);
1601 if (bNeedRepartition && DOMAINDECOMP(cr))
1603 dd_collect_state(cr->dd, state, state_global);
1607 /* Replica exchange */
1608 bExchanged = FALSE;
1609 if (bDoReplEx)
1611 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1614 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1616 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1, state_global, *top_global, ir,
1617 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
1618 nrnb, wcycle, FALSE);
1619 shouldCheckNumberOfBondedInteractions = true;
1620 upd.setNumAtoms(state->natoms);
1623 bFirstStep = FALSE;
1624 bInitStep = FALSE;
1626 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1627 /* With all integrators, except VV, we need to retain the pressure
1628 * at the current step for coupling at the next step.
1630 if ((state->flags & (1U << estPRES_PREV))
1631 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1633 /* Store the pressure in t_state for pressure coupling
1634 * at the next MD step.
1636 copy_mat(pres, state->pres_prev);
1639 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1641 if ((membed != nullptr) && (!bLastStep))
1643 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1646 cycles = wallcycle_stop(wcycle, ewcSTEP);
1647 if (DOMAINDECOMP(cr) && wcycle)
1649 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1652 /* increase the MD step number */
1653 step++;
1654 step_rel++;
1656 resetHandler->resetCounters(step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb,
1657 fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1659 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1660 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1662 /* End of main MD loop */
1664 /* Closing TNG files can include compressing data. Therefore it is good to do that
1665 * before stopping the time measurements. */
1666 mdoutf_tng_close(outf);
1668 /* Stop measuring walltime */
1669 walltime_accounting_end_time(walltime_accounting);
1671 if (!thisRankHasDuty(cr, DUTY_PME))
1673 /* Tell the PME only node to finish */
1674 gmx_pme_send_finish(cr);
1677 if (MASTER(cr))
1679 if (ir->nstcalcenergy > 0)
1681 gmx::EnergyOutput::printAnnealingTemperatures(fplog, groups, &(ir->opts));
1682 energyOutput.printAverages(fplog, groups);
1685 done_mdoutf(outf);
1687 if (bPMETune)
1689 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1692 done_shellfc(fplog, shellfc, step_rel);
1694 if (useReplicaExchange && MASTER(cr))
1696 print_replica_exchange_statistics(fplog, repl_ex);
1699 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1701 global_stat_destroy(gstat);