Fix cycle counters for "comm.coord" and "Wait + Comm. F" to support GPU halo exchange...
[gromacs.git] / src / gromacs / mdrun / md.cpp
blob5f5f863e816b40c8f93b9898ee1ccc67dde07a62
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/listed_forces.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/coupling.h"
82 #include "gromacs/mdlib/ebin.h"
83 #include "gromacs/mdlib/enerdata_utils.h"
84 #include "gromacs/mdlib/energyoutput.h"
85 #include "gromacs/mdlib/expanded.h"
86 #include "gromacs/mdlib/force.h"
87 #include "gromacs/mdlib/force_flags.h"
88 #include "gromacs/mdlib/forcerec.h"
89 #include "gromacs/mdlib/md_support.h"
90 #include "gromacs/mdlib/mdatoms.h"
91 #include "gromacs/mdlib/mdoutf.h"
92 #include "gromacs/mdlib/membed.h"
93 #include "gromacs/mdlib/resethandler.h"
94 #include "gromacs/mdlib/sighandler.h"
95 #include "gromacs/mdlib/simulationsignal.h"
96 #include "gromacs/mdlib/stat.h"
97 #include "gromacs/mdlib/stophandler.h"
98 #include "gromacs/mdlib/tgroup.h"
99 #include "gromacs/mdlib/trajectory_writing.h"
100 #include "gromacs/mdlib/update.h"
101 #include "gromacs/mdlib/update_constrain_gpu.h"
102 #include "gromacs/mdlib/vcm.h"
103 #include "gromacs/mdlib/vsite.h"
104 #include "gromacs/mdrunutility/handlerestart.h"
105 #include "gromacs/mdrunutility/multisim.h"
106 #include "gromacs/mdrunutility/printtime.h"
107 #include "gromacs/mdtypes/awh_history.h"
108 #include "gromacs/mdtypes/awh_params.h"
109 #include "gromacs/mdtypes/commrec.h"
110 #include "gromacs/mdtypes/df_history.h"
111 #include "gromacs/mdtypes/energyhistory.h"
112 #include "gromacs/mdtypes/fcdata.h"
113 #include "gromacs/mdtypes/forcerec.h"
114 #include "gromacs/mdtypes/group.h"
115 #include "gromacs/mdtypes/inputrec.h"
116 #include "gromacs/mdtypes/interaction_const.h"
117 #include "gromacs/mdtypes/md_enums.h"
118 #include "gromacs/mdtypes/mdatom.h"
119 #include "gromacs/mdtypes/mdrunoptions.h"
120 #include "gromacs/mdtypes/observableshistory.h"
121 #include "gromacs/mdtypes/pullhistory.h"
122 #include "gromacs/mdtypes/simulation_workload.h"
123 #include "gromacs/mdtypes/state.h"
124 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
125 #include "gromacs/modularsimulator/energyelement.h"
126 #include "gromacs/nbnxm/gpu_data_mgmt.h"
127 #include "gromacs/nbnxm/nbnxm.h"
128 #include "gromacs/pbcutil/pbc.h"
129 #include "gromacs/pulling/output.h"
130 #include "gromacs/pulling/pull.h"
131 #include "gromacs/swap/swapcoords.h"
132 #include "gromacs/timing/wallcycle.h"
133 #include "gromacs/timing/walltime_accounting.h"
134 #include "gromacs/topology/atoms.h"
135 #include "gromacs/topology/idef.h"
136 #include "gromacs/topology/mtop_util.h"
137 #include "gromacs/topology/topology.h"
138 #include "gromacs/trajectory/trajectoryframe.h"
139 #include "gromacs/utility/basedefinitions.h"
140 #include "gromacs/utility/cstringutil.h"
141 #include "gromacs/utility/fatalerror.h"
142 #include "gromacs/utility/logger.h"
143 #include "gromacs/utility/real.h"
144 #include "gromacs/utility/smalloc.h"
146 #include "legacysimulator.h"
147 #include "replicaexchange.h"
148 #include "shellfc.h"
150 #if GMX_FAHCORE
151 # include "corewrap.h"
152 #endif
154 using gmx::SimulationSignaller;
156 void gmx::LegacySimulator::do_md()
158 // TODO Historically, the EM and MD "integrators" used different
159 // names for the t_inputrec *parameter, but these must have the
160 // same name, now that it's a member of a struct. We use this ir
161 // alias to avoid a large ripple of nearly useless changes.
162 // t_inputrec is being replaced by IMdpOptionsProvider, so this
163 // will go away eventually.
164 t_inputrec* ir = inputrec;
165 int64_t step, step_rel;
166 double t, t0 = ir->init_t, lam0[efptNR];
167 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
168 gmx_bool bNS, bNStList, bStopCM, bFirstStep, bInitStep, bLastStep = FALSE;
169 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
170 gmx_bool do_ene, do_log, do_verbose;
171 gmx_bool bMasterState;
172 unsigned int force_flags;
173 tensor force_vir = { { 0 } }, shake_vir = { { 0 } }, total_vir = { { 0 } }, 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 const t_fcdata& fcdata = fr->listedForces->fcdata();
263 bool simulationsShareState = false;
264 int nstSignalComm = nstglobalcomm;
266 // TODO This implementation of ensemble orientation restraints is nasty because
267 // a user can't just do multi-sim with single-sim orientation restraints.
268 bool usingEnsembleRestraints =
269 (fcdata.disres->nsystems > 1) || ((ms != nullptr) && (fcdata.orires->nr != 0));
270 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
272 // Replica exchange, ensemble restraints and AWH need all
273 // simulations to remain synchronized, so they need
274 // checkpoints and stop conditions to act on the same step, so
275 // the propagation of such signals must take place between
276 // simulations, not just within simulations.
277 // TODO: Make algorithm initializers set these flags.
278 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
280 if (simulationsShareState)
282 // Inter-simulation signal communication does not need to happen
283 // often, so we use a minimum of 200 steps to reduce overhead.
284 const int c_minimumInterSimulationSignallingInterval = 200;
285 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
286 * nstglobalcomm;
290 if (startingBehavior != StartingBehavior::RestartWithAppending)
292 pleaseCiteCouplingAlgorithms(fplog, *ir);
294 gmx_mdoutf* outf =
295 init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, ir,
296 top_global, oenv, wcycle, startingBehavior, simulationsShareState, ms);
297 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
298 mdoutf_get_fp_dhdl(outf), false, startingBehavior, mdModulesNotifier);
300 gstat = global_stat_init(ir);
302 /* Check for polarizable models and flexible constraints */
303 shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
304 ir->nstcalcenergy, DOMAINDECOMP(cr));
307 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
308 if ((io > 2000) && MASTER(cr))
310 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
314 // Local state only becomes valid now.
315 std::unique_ptr<t_state> stateInstance;
316 t_state* state;
319 gmx_localtop_t top(top_global->ffparams);
321 auto mdatoms = mdAtoms->mdatoms();
323 std::unique_ptr<UpdateConstrainGpu> integrator;
325 if (DOMAINDECOMP(cr))
327 stateInstance = std::make_unique<t_state>();
328 state = stateInstance.get();
329 dd_init_local_state(cr->dd, state_global, state);
331 /* Distribute the charge groups over the nodes from the master node */
332 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
333 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
334 nrnb, nullptr, FALSE);
335 shouldCheckNumberOfBondedInteractions = true;
336 upd.setNumAtoms(state->natoms);
338 else
340 state_change_natoms(state_global, state_global->natoms);
341 /* Copy the pointer to the global state */
342 state = state_global;
344 /* Generate and initialize new topology */
345 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
347 upd.setNumAtoms(state->natoms);
350 const auto& simulationWork = runScheduleWork->simulationWork;
351 const bool useGpuForPme = simulationWork.useGpuPme;
352 const bool useGpuForNonbonded = simulationWork.useGpuNonbonded;
353 const bool useGpuForBufferOps = simulationWork.useGpuBufferOps;
354 const bool useGpuForUpdate = simulationWork.useGpuUpdate;
356 StatePropagatorDataGpu* stateGpu = fr->stateGpu;
358 // TODO: the assertions below should be handled by UpdateConstraintsBuilder.
359 if (useGpuForUpdate)
361 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
362 || constr->numConstraintsTotal() == 0,
363 "Constraints in domain decomposition are only supported with update "
364 "groups if using GPU update.\n");
365 GMX_RELEASE_ASSERT(ir->eConstrAlg != econtSHAKE || constr == nullptr
366 || constr->numConstraintsTotal() == 0,
367 "SHAKE is not supported with GPU update.");
368 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
369 "Either PME or short-ranged non-bonded interaction tasks must run on "
370 "the GPU to use GPU update.\n");
371 GMX_RELEASE_ASSERT(ir->eI == eiMD,
372 "Only the md integrator is supported with the GPU update.\n");
373 GMX_RELEASE_ASSERT(
374 ir->etc != etcNOSEHOOVER,
375 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
376 GMX_RELEASE_ASSERT(ir->epc == epcNO || ir->epc == epcPARRINELLORAHMAN || ir->epc == epcBERENDSEN,
377 "Only Parrinello-Rahman and Berendsen pressure coupling are supported "
378 "with the GPU update.\n");
379 GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
380 "Virtual sites are not supported with the GPU update.\n");
381 GMX_RELEASE_ASSERT(ed == nullptr,
382 "Essential dynamics is not supported with the GPU update.\n");
383 GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(ir->pull),
384 "Constraints pulling is not supported with the GPU update.\n");
385 GMX_RELEASE_ASSERT(fcdata.orires->nr == 0,
386 "Orientation restraints are not supported with the GPU update.\n");
387 GMX_RELEASE_ASSERT(
388 ir->efep == efepNO
389 || (!haveFreeEnergyType(*ir, efptBONDED) && !haveFreeEnergyType(*ir, efptMASS)),
390 "Free energy perturbation of masses and constraints are not supported with the GPU "
391 "update.");
393 if (constr != nullptr && constr->numConstraintsTotal() > 0)
395 GMX_LOG(mdlog.info)
396 .asParagraph()
397 .appendText("Updating coordinates and applying constraints on the GPU.");
399 else
401 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
403 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
404 "Device stream manager should be initialized in order to use GPU "
405 "update-constraints.");
406 GMX_RELEASE_ASSERT(
407 fr->deviceStreamManager->streamIsValid(gmx::DeviceStreamType::UpdateAndConstraints),
408 "Update stream should be initialized in order to use GPU "
409 "update-constraints.");
410 integrator = std::make_unique<UpdateConstrainGpu>(
411 *ir, *top_global, fr->deviceStreamManager->context(),
412 fr->deviceStreamManager->stream(gmx::DeviceStreamType::UpdateAndConstraints),
413 stateGpu->xUpdatedOnDevice());
415 integrator->setPbc(PbcType::Xyz, state->box);
418 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
420 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
422 if ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
424 changePinningPolicy(&f, PinningPolicy::PinnedIfSupported);
426 if (useGpuForUpdate)
428 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
431 // NOTE: The global state is no longer used at this point.
432 // But state_global is still used as temporary storage space for writing
433 // the global state to file and potentially for replica exchange.
434 // (Global topology should persist.)
436 update_mdatoms(mdatoms, state->lambda[efptMASS]);
438 if (ir->bExpanded)
440 /* Check nstexpanded here, because the grompp check was broken */
441 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
443 gmx_fatal(FARGS,
444 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
446 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
449 if (MASTER(cr))
451 EnergyElement::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
454 preparePrevStepPullCom(ir, pull_work, mdatoms->massT, state, state_global, cr,
455 startingBehavior != StartingBehavior::NewSimulation);
457 // TODO: Remove this by converting AWH into a ForceProvider
458 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms,
459 startingBehavior != StartingBehavior::NewSimulation,
460 shellfc != nullptr, opt2fn("-awh", nfile, fnm), pull_work);
462 if (useReplicaExchange && MASTER(cr))
464 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir, replExParams);
466 /* PME tuning is only supported in the Verlet scheme, with PME for
467 * Coulomb. It is not supported with only LJ PME. */
468 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
469 && ir->cutoff_scheme != ecutsGROUP);
471 pme_load_balancing_t* pme_loadbal = nullptr;
472 if (bPMETune)
474 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata,
475 fr->nbv->useGpu());
478 if (!ir->bContinuation)
480 if (state->flags & (1U << estV))
482 auto v = makeArrayRef(state->v);
483 /* Set the velocities of vsites, shells and frozen atoms to zero */
484 for (i = 0; i < mdatoms->homenr; i++)
486 if (mdatoms->ptype[i] == eptVSite || mdatoms->ptype[i] == eptShell)
488 clear_rvec(v[i]);
490 else if (mdatoms->cFREEZE)
492 for (m = 0; m < DIM; m++)
494 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
496 v[i][m] = 0;
503 if (constr)
505 /* Constrain the initial coordinates and velocities */
506 do_constrain_first(fplog, constr, ir, mdatoms->nr, mdatoms->homenr,
507 state->x.arrayRefWithPadding(), state->v.arrayRefWithPadding(),
508 state->box, state->lambda[efptBONDED]);
510 if (vsite)
512 /* Construct the virtual sites for the initial configuration */
513 vsite->construct(state->x, ir->delta_t, {}, state->box);
517 if (ir->efep != efepNO)
519 /* Set free energy calculation frequency as the greatest common
520 * denominator of nstdhdl and repl_ex_nst. */
521 nstfep = ir->fepvals->nstdhdl;
522 if (ir->bExpanded)
524 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
526 if (useReplicaExchange)
528 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
532 /* Be REALLY careful about what flags you set here. You CANNOT assume
533 * this is the first step, since we might be restarting from a checkpoint,
534 * and in that case we should not do any modifications to the state.
536 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
538 // When restarting from a checkpoint, it can be appropriate to
539 // initialize ekind from quantities in the checkpoint. Otherwise,
540 // compute_globals must initialize ekind before the simulation
541 // starts/restarts. However, only the master rank knows what was
542 // found in the checkpoint file, so we have to communicate in
543 // order to coordinate the restart.
545 // TODO Consider removing this communication if/when checkpoint
546 // reading directly follows .tpr reading, because all ranks can
547 // agree on hasReadEkinState at that time.
548 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
549 if (PAR(cr))
551 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr->mpi_comm_mygroup);
553 if (hasReadEkinState)
555 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
558 unsigned int cglo_flags =
559 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
560 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
562 bSumEkinhOld = FALSE;
564 t_vcm vcm(top_global->groups, *ir);
565 reportComRemovalInfo(fplog, vcm);
567 /* To minimize communication, compute_globals computes the COM velocity
568 * and the kinetic energy for the velocities without COM motion removed.
569 * Thus to get the kinetic energy without the COM contribution, we need
570 * to call compute_globals twice.
572 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
574 unsigned int cglo_flags_iteration = cglo_flags;
575 if (bStopCM && cgloIteration == 0)
577 cglo_flags_iteration |= CGLO_STOPCM;
578 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
580 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
581 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, nullptr,
582 enerd, force_vir, shake_vir, total_vir, pres, constr, &nullSignaller,
583 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
584 cglo_flags_iteration
585 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
586 : 0));
587 if (cglo_flags_iteration & CGLO_STOPCM)
589 /* At initialization, do not pass x with acceleration-correction mode
590 * to avoid (incorrect) correction of the initial coordinates.
592 auto x = (vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION) ? ArrayRef<RVec>()
593 : makeArrayRef(state->x);
594 process_and_stopcm_grp(fplog, &vcm, *mdatoms, x, makeArrayRef(state->v));
595 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
598 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
599 makeConstArrayRef(state->x), state->box,
600 &shouldCheckNumberOfBondedInteractions);
601 if (ir->eI == eiVVAK)
603 /* a second call to get the half step temperature initialized as well */
604 /* we do the same call as above, but turn the pressure off -- internally to
605 compute_globals, this is recognized as a velocity verlet half-step
606 kinetic energy calculation. This minimized excess variables, but
607 perhaps loses some logic?*/
609 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
610 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, nullptr,
611 enerd, force_vir, shake_vir, total_vir, pres, constr, &nullSignaller,
612 state->box, nullptr, &bSumEkinhOld, cglo_flags & ~CGLO_PRESSURE);
615 /* Calculate the initial half step temperature, and save the ekinh_old */
616 if (startingBehavior == StartingBehavior::NewSimulation)
618 for (i = 0; (i < ir->opts.ngtc); i++)
620 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
624 /* need to make an initiation call to get the Trotter variables set, as well as other constants
625 for non-trotter temperature control */
626 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
628 if (MASTER(cr))
630 if (!ir->bContinuation)
632 if (constr && ir->eConstrAlg == econtLINCS)
634 fprintf(fplog, "RMS relative constraint deviation after constraining: %.2e\n",
635 constr->rmsd());
637 if (EI_STATE_VELOCITY(ir->eI))
639 real temp = enerd->term[F_TEMP];
640 if (ir->eI != eiVV)
642 /* Result of Ekin averaged over velocities of -half
643 * and +half step, while we only have -half step here.
645 temp *= 2;
647 fprintf(fplog, "Initial temperature: %g K\n", temp);
651 char tbuf[20];
652 fprintf(stderr, "starting mdrun '%s'\n", *(top_global->name));
653 if (ir->nsteps >= 0)
655 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
657 else
659 sprintf(tbuf, "%s", "infinite");
661 if (ir->init_step > 0)
663 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
664 gmx_step_str(ir->init_step + ir->nsteps, sbuf), tbuf,
665 gmx_step_str(ir->init_step, sbuf2), ir->init_step * ir->delta_t);
667 else
669 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
671 fprintf(fplog, "\n");
674 walltime_accounting_start_time(walltime_accounting);
675 wallcycle_start(wcycle, ewcRUN);
676 print_start(fplog, cr, walltime_accounting, "mdrun");
678 #if GMX_FAHCORE
679 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
680 int chkpt_ret = fcCheckPointParallel(cr->nodeid, NULL, 0);
681 if (chkpt_ret == 0)
683 gmx_fatal(3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0);
685 #endif
687 /***********************************************************
689 * Loop over MD steps
691 ************************************************************/
693 bFirstStep = TRUE;
694 /* Skip the first Nose-Hoover integration when we get the state from tpx */
695 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
696 bSumEkinhOld = FALSE;
697 bExchanged = FALSE;
698 bNeedRepartition = FALSE;
700 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
701 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
702 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
703 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
705 auto checkpointHandler = std::make_unique<CheckpointHandler>(
706 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]), simulationsShareState,
707 ir->nstlist == 0, MASTER(cr), mdrunOptions.writeConfout,
708 mdrunOptions.checkpointOptions.period);
710 const bool resetCountersIsLocal = true;
711 auto resetHandler = std::make_unique<ResetHandler>(
712 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
713 !resetCountersIsLocal, ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
714 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
716 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
718 step = ir->init_step;
719 step_rel = 0;
721 // TODO extract this to new multi-simulation module
722 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
724 if (!multisim_int_all_are_equal(ms, ir->nsteps))
726 GMX_LOG(mdlog.warning)
727 .appendText(
728 "Note: The number of steps is not consistent across multi "
729 "simulations,\n"
730 "but we are proceeding anyway!");
732 if (!multisim_int_all_are_equal(ms, ir->init_step))
734 if (simulationsShareState)
736 if (MASTER(cr))
738 gmx_fatal(FARGS,
739 "The initial step is not consistent across multi simulations which "
740 "share the state");
742 gmx_barrier(cr->mpi_comm_mygroup);
744 else
746 GMX_LOG(mdlog.warning)
747 .appendText(
748 "Note: The initial step is not consistent across multi "
749 "simulations,\n"
750 "but we are proceeding anyway!");
755 /* and stop now if we should */
756 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
757 while (!bLastStep)
760 /* Determine if this is a neighbor search step */
761 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
763 if (bPMETune && bNStList)
765 // This has to be here because PME load balancing is called so early.
766 // TODO: Move to after all booleans are defined.
767 if (useGpuForUpdate && !bFirstStep)
769 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
770 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
772 /* PME grid + cut-off optimization with GPUs or PME nodes */
773 pme_loadbal_do(pme_loadbal, cr, (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
774 fplog, mdlog, *ir, fr, state->box, state->x, wcycle, step, step_rel,
775 &bPMETunePrinting, simulationWork.useGpuPmePpCommunication);
778 wallcycle_start(wcycle, ewcSTEP);
780 bLastStep = (step_rel == ir->nsteps);
781 t = t0 + step * ir->delta_t;
783 // TODO Refactor this, so that nstfep does not need a default value of zero
784 if (ir->efep != efepNO || ir->bSimTemp)
786 /* find and set the current lambdas */
787 setCurrentLambdasLocal(step, ir->fepvals, lam0, state->lambda, state->fep_state);
789 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
790 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
791 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
792 && (!bFirstStep));
795 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
796 && do_per_step(step, replExParams.exchangeInterval));
798 if (doSimulatedAnnealing)
800 update_annealing_target_temp(ir, t, &upd);
803 /* Stop Center of Mass motion */
804 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
806 /* Determine whether or not to do Neighbour Searching */
807 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
809 /* Note that the stopHandler will cause termination at nstglobalcomm
810 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
811 * nstpcouple steps, we have computed the half-step kinetic energy
812 * of the previous step and can always output energies at the last step.
814 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
816 /* do_log triggers energy and virial calculation. Because this leads
817 * to different code paths, forces can be different. Thus for exact
818 * continuation we should avoid extra log output.
819 * Note that the || bLastStep can result in non-exact continuation
820 * beyond the last step. But we don't consider that to be an issue.
822 do_log = (do_per_step(step, ir->nstlog)
823 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
824 do_verbose = mdrunOptions.verbose
825 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
827 if (useGpuForUpdate && !bFirstStep && bNS)
829 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
830 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
831 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
832 // Copy coordinate from the GPU when needed at the search step.
833 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
834 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
835 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
836 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
839 if (bNS && !(bFirstStep && ir->bContinuation))
841 bMasterState = FALSE;
842 /* Correct the new box if it is too skewed */
843 if (inputrecDynamicBox(ir))
845 if (correct_box(fplog, step, state->box))
847 bMasterState = TRUE;
848 // If update is offloaded, it should be informed about the box size change
849 if (useGpuForUpdate)
851 integrator->setPbc(PbcType::Xyz, state->box);
855 if (DOMAINDECOMP(cr) && bMasterState)
857 dd_collect_state(cr->dd, state, state_global);
860 if (DOMAINDECOMP(cr))
862 /* Repartition the domain decomposition */
863 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
864 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
865 fr, vsite, constr, nrnb, wcycle, do_verbose && !bPMETunePrinting);
866 shouldCheckNumberOfBondedInteractions = true;
867 upd.setNumAtoms(state->natoms);
869 // Allocate or re-size GPU halo exchange object, if necessary
870 if (havePPDomainDecomposition(cr) && simulationWork.useGpuHaloExchange
871 && useGpuForNonbonded && is1D(*cr->dd))
873 GMX_RELEASE_ASSERT(fr->deviceStreamManager != nullptr,
874 "GPU device manager has to be initialized to use GPU "
875 "version of halo exchange.");
876 // TODO remove need to pass local stream into GPU halo exchange - Issue #3093
877 constructGpuHaloExchange(mdlog, *cr, *fr->deviceStreamManager, wcycle);
882 if (MASTER(cr) && do_log)
884 gmx::EnergyOutput::printHeader(fplog, step,
885 t); /* can we improve the information printed here? */
888 if (ir->efep != efepNO)
890 update_mdatoms(mdatoms, state->lambda[efptMASS]);
893 if (bExchanged)
896 /* We need the kinetic energy at minus the half step for determining
897 * the full step kinetic energy and possibly for T-coupling.*/
898 /* This may not be quite working correctly yet . . . . */
899 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
900 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, wcycle,
901 enerd, nullptr, nullptr, nullptr, nullptr, constr, &nullSignaller,
902 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
903 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
904 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
905 &top, makeConstArrayRef(state->x), state->box,
906 &shouldCheckNumberOfBondedInteractions);
908 clear_mat(force_vir);
910 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
912 /* Determine the energy and pressure:
913 * at nstcalcenergy steps and at energy output steps (set below).
915 if (EI_VV(ir->eI) && (!bInitStep))
917 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
918 bCalcVir = bCalcEnerStep
919 || (ir->epc != epcNO
920 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
922 else
924 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
925 bCalcVir = bCalcEnerStep || (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
927 bCalcEner = bCalcEnerStep;
929 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
931 if (do_ene || do_log || bDoReplEx)
933 bCalcVir = TRUE;
934 bCalcEner = TRUE;
937 /* Do we need global communication ? */
938 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
939 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
941 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
942 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
943 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
945 if (shellfc)
947 /* Now is the time to relax the shells */
948 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
949 imdSession, pull_work, bNS, force_flags, &top, constr, enerd,
950 state->natoms, state->x.arrayRefWithPadding(),
951 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
952 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, shellfc,
953 fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
955 else
957 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
958 is updated (or the AWH update will be performed twice for one step when continuing).
959 It would be best to call this update function from do_md_trajectory_writing but that
960 would occur after do_force. One would have to divide the update_awh function into one
961 function applying the AWH force and one doing the AWH bias update. The update AWH
962 bias function could then be called after do_md_trajectory_writing (then containing
963 update_awh_history). The checkpointing will in the future probably moved to the start
964 of the md loop which will rid of this issue. */
965 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
967 awh->updateHistory(state_global->awhHistory.get());
970 /* The coordinates (x) are shifted (to get whole molecules)
971 * in do_force.
972 * This is parallellized as well, and does communication too.
973 * Check comments in sim_util.c
975 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation, imdSession, pull_work, step,
976 nrnb, wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
977 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, state->lambda, fr,
978 runScheduleWork, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
979 (bNS ? GMX_FORCE_NS : 0) | force_flags, ddBalanceRegionHandler);
982 // VV integrators do not need the following velocity half step
983 // if it is the first step after starting from a checkpoint.
984 // That is, the half step is needed on all other steps, and
985 // also the first step when starting from a .tpr file.
986 if (EI_VV(ir->eI) && (!bFirstStep || startingBehavior == StartingBehavior::NewSimulation))
987 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
989 rvec* vbuf = nullptr;
991 wallcycle_start(wcycle, ewcUPDATE);
992 if (ir->eI == eiVV && bInitStep)
994 /* if using velocity verlet with full time step Ekin,
995 * take the first half step only to compute the
996 * virial for the first step. From there,
997 * revert back to the initial coordinates
998 * so that the input is actually the initial step.
1000 snew(vbuf, state->natoms);
1001 copy_rvecn(state->v.rvec_array(), vbuf, 0,
1002 state->natoms); /* should make this better for parallelizing? */
1004 else
1006 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1007 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1008 trotter_seq, ettTSEQ1);
1011 upd.update_coords(*ir, step, mdatoms, state, f.arrayRefWithPadding(), fcdata, ekind, M,
1012 etrtVELOCITY1, cr, constr != nullptr);
1014 wallcycle_stop(wcycle, ewcUPDATE);
1015 constrain_velocities(constr, do_log, do_ene, step, state, nullptr, bCalcVir, shake_vir);
1016 wallcycle_start(wcycle, ewcUPDATE);
1017 /* if VV, compute the pressure and constraints */
1018 /* For VV2, we strictly only need this if using pressure
1019 * control, but we really would like to have accurate pressures
1020 * printed out.
1021 * Think about ways around this in the future?
1022 * For now, keep this choice in comments.
1024 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1025 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1026 bPres = TRUE;
1027 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1028 if (bCalcEner && ir->eI == eiVVAK)
1030 bSumEkinhOld = TRUE;
1032 /* for vv, the first half of the integration actually corresponds to the previous step.
1033 So we need information from the last step in the first half of the integration */
1034 if (bGStat || do_per_step(step - 1, nstglobalcomm))
1036 wallcycle_stop(wcycle, ewcUPDATE);
1037 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1038 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, wcycle,
1039 enerd, force_vir, shake_vir, total_vir, pres, constr, &nullSignaller,
1040 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1041 (bGStat ? CGLO_GSTAT : 0) | (bCalcEner ? CGLO_ENERGY : 0)
1042 | (bTemp ? CGLO_TEMPERATURE : 0) | (bPres ? CGLO_PRESSURE : 0)
1043 | (bPres ? CGLO_CONSTRAINT : 0) | (bStopCM ? CGLO_STOPCM : 0)
1044 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1045 : 0)
1046 | CGLO_SCALEEKIN);
1047 /* explanation of above:
1048 a) We compute Ekin at the full time step
1049 if 1) we are using the AveVel Ekin, and it's not the
1050 initial step, or 2) if we are using AveEkin, but need the full
1051 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1052 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1053 EkinAveVel because it's needed for the pressure */
1054 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1055 top_global, &top, makeConstArrayRef(state->x),
1056 state->box, &shouldCheckNumberOfBondedInteractions);
1057 if (bStopCM)
1059 process_and_stopcm_grp(fplog, &vcm, *mdatoms, makeArrayRef(state->x),
1060 makeArrayRef(state->v));
1061 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1063 wallcycle_start(wcycle, ewcUPDATE);
1065 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1066 if (!bInitStep)
1068 if (bTrotter)
1070 m_add(force_vir, shake_vir,
1071 total_vir); /* we need the un-dispersion corrected total vir here */
1072 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1073 trotter_seq, ettTSEQ2);
1075 /* TODO This is only needed when we're about to write
1076 * a checkpoint, because we use it after the restart
1077 * (in a kludge?). But what should we be doing if
1078 * the startingBehavior is NewSimulation or bInitStep are true? */
1079 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1081 copy_mat(shake_vir, state->svir_prev);
1082 copy_mat(force_vir, state->fvir_prev);
1084 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1086 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1087 enerd->term[F_TEMP] =
1088 sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1089 enerd->term[F_EKIN] = trace(ekind->ekin);
1092 else if (bExchanged)
1094 wallcycle_stop(wcycle, ewcUPDATE);
1095 /* We need the kinetic energy at minus the half step for determining
1096 * the full step kinetic energy and possibly for T-coupling.*/
1097 /* This may not be quite working correctly yet . . . . */
1098 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1099 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm, wcycle,
1100 enerd, nullptr, nullptr, nullptr, nullptr, constr, &nullSignaller,
1101 state->box, nullptr, &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(constr, do_log, do_ene, step, state, nullptr, false, nullptr);
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 upd.update_coords(*ir, step, mdatoms, state, f.arrayRefWithPadding(), fcdata, ekind, M,
1260 etrtVELOCITY2, cr, constr != nullptr);
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 upd.update_coords(*ir, step, mdatoms, state, f.arrayRefWithPadding(), fcdata, ekind, M,
1326 etrtPOSITION, cr, constr != nullptr);
1328 wallcycle_stop(wcycle, ewcUPDATE);
1330 constrain_coordinates(constr, do_log, do_ene, step, state,
1331 upd.xp()->arrayRefWithPadding(), &dvdl_constr, bCalcVir, shake_vir);
1333 upd.update_sd_second_half(*ir, step, &dvdl_constr, mdatoms, state, cr, nrnb, wcycle,
1334 constr, do_log, do_ene);
1335 upd.finish_update(*ir, mdatoms, state, wcycle, constr != nullptr);
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, mdatoms, nrnb, &vcm, wcycle, enerd,
1349 force_vir, shake_vir, total_vir, pres, constr, &nullSignaller, lastbox,
1350 nullptr, &bSumEkinhOld, (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE);
1351 wallcycle_start(wcycle, ewcUPDATE);
1352 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1353 /* now we know the scaling, we can compute the positions again */
1354 std::copy(cbuf.begin(), cbuf.end(), state->x.begin());
1356 upd.update_coords(*ir, step, mdatoms, state, f.arrayRefWithPadding(), fcdata, ekind, M,
1357 etrtPOSITION, cr, constr != nullptr);
1358 wallcycle_stop(wcycle, ewcUPDATE);
1360 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1361 /* are the small terms in the shake_vir here due
1362 * to numerical errors, or are they important
1363 * physically? I'm thinking they are just errors, but not completely sure.
1364 * For now, will call without actually constraining, constr=NULL*/
1365 upd.finish_update(*ir, mdatoms, state, wcycle, false);
1367 if (EI_VV(ir->eI))
1369 /* this factor or 2 correction is necessary
1370 because half of the constraint force is removed
1371 in the vv step, so we have to double it. See
1372 the Issue #1255. It is not yet clear
1373 if the factor of 2 is exact, or just a very
1374 good approximation, and this will be
1375 investigated. The next step is to see if this
1376 can be done adding a dhdl contribution from the
1377 rattle step, but this is somewhat more
1378 complicated with the current code. Will be
1379 investigated, hopefully for 4.6.3. However,
1380 this current solution is much better than
1381 having it completely wrong.
1383 enerd->term[F_DVDL_CONSTR] += 2 * dvdl_constr;
1385 else
1387 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1390 if (vsite != nullptr)
1392 wallcycle_start(wcycle, ewcVSITECONSTR);
1393 vsite->construct(state->x, ir->delta_t, state->v, state->box);
1394 wallcycle_stop(wcycle, ewcVSITECONSTR);
1397 /* ############## IF NOT VV, Calculate globals HERE ############ */
1398 /* With Leap-Frog we can skip compute_globals at
1399 * non-communication steps, but we need to calculate
1400 * the kinetic energy one step before communication.
1403 // Organize to do inter-simulation signalling on steps if
1404 // and when algorithms require it.
1405 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1407 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1409 // Copy coordinates when needed to stop the CM motion.
1410 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1412 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1413 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1415 // Since we're already communicating at this step, we
1416 // can propagate intra-simulation signals. Note that
1417 // check_nstglobalcomm has the responsibility for
1418 // choosing the value of nstglobalcomm that is one way
1419 // bGStat becomes true, so we can't get into a
1420 // situation where e.g. checkpointing can't be
1421 // signalled.
1422 bool doIntraSimSignal = true;
1423 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1425 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
1426 makeConstArrayRef(state->v), state->box, mdatoms, nrnb, &vcm,
1427 wcycle, enerd, force_vir, shake_vir, total_vir, pres, constr,
1428 &signaller, lastbox, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1429 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1430 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1431 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1432 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1433 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1434 : 0));
1435 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1436 top_global, &top, makeConstArrayRef(state->x),
1437 state->box, &shouldCheckNumberOfBondedInteractions);
1438 if (!EI_VV(ir->eI) && bStopCM)
1440 process_and_stopcm_grp(fplog, &vcm, *mdatoms, makeArrayRef(state->x),
1441 makeArrayRef(state->v));
1442 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1444 // TODO: The special case of removing CM motion should be dealt more gracefully
1445 if (useGpuForUpdate)
1447 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1448 // Here we block until the H2D copy completes because event sync with the
1449 // force kernels that use the coordinates on the next steps is not implemented
1450 // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1451 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1452 // If the COM removal changed the velocities on the CPU, this has to be accounted for.
1453 if (vcm.mode != ecmNO)
1455 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1462 /* ############# END CALC EKIN AND PRESSURE ################# */
1464 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1465 the virial that should probably be addressed eventually. state->veta has better properies,
1466 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1467 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1469 if (ir->efep != efepNO && !EI_VV(ir->eI))
1471 /* Sum up the foreign energy and dhdl terms for md and sd.
1472 Currently done every step so that dhdl is correct in the .edr */
1473 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1476 update_pcouple_after_coordinates(fplog, step, ir, mdatoms, pres, force_vir, shake_vir,
1477 pressureCouplingMu, state, nrnb, upd.deform(), !useGpuForUpdate);
1479 const bool doBerendsenPressureCoupling =
1480 (inputrec->epc == epcBERENDSEN && do_per_step(step, inputrec->nstpcouple));
1481 if (useGpuForUpdate && (doBerendsenPressureCoupling || doParrinelloRahman))
1483 integrator->scaleCoordinates(pressureCouplingMu);
1484 integrator->setPbc(PbcType::Xyz, state->box);
1487 /* ################# END UPDATE STEP 2 ################# */
1488 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1490 /* The coordinates (x) were unshifted in update */
1491 if (!bGStat)
1493 /* We will not sum ekinh_old,
1494 * so signal that we still have to do it.
1496 bSumEkinhOld = TRUE;
1499 if (bCalcEner)
1501 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1503 /* use the directly determined last velocity, not actually the averaged half steps */
1504 if (bTrotter && ir->eI == eiVV)
1506 enerd->term[F_EKIN] = last_ekin;
1508 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1510 if (integratorHasConservedEnergyQuantity(ir))
1512 if (EI_VV(ir->eI))
1514 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1516 else
1518 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1521 /* ######### END PREPARING EDR OUTPUT ########### */
1524 /* Output stuff */
1525 if (MASTER(cr))
1527 if (fplog && do_log && bDoExpanded)
1529 /* only needed if doing expanded ensemble */
1530 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals,
1531 ir->bSimTemp ? ir->simtempvals : nullptr,
1532 state_global->dfhist, state->fep_state, ir->nstlog, step);
1534 if (bCalcEner)
1536 energyOutput.addDataAtEnergyStep(bDoDHDL, bCalcEnerStep, t, mdatoms->tmass, enerd, state,
1537 ir->fepvals, ir->expandedvals, lastbox, shake_vir,
1538 force_vir, total_vir, pres, ekind, mu_tot, constr);
1540 else
1542 energyOutput.recordNonEnergyStep();
1545 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1546 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1548 if (doSimulatedAnnealing)
1550 gmx::EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups,
1551 &(ir->opts));
1553 if (do_log || do_ene || do_dr || do_or)
1555 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1556 do_log ? fplog : nullptr, step, t,
1557 &fr->listedForces->fcdata(), awh.get());
1560 if (ir->bPull)
1562 pull_print_output(pull_work, step, t);
1565 if (do_per_step(step, ir->nstlog))
1567 if (fflush(fplog) != 0)
1569 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1573 if (bDoExpanded)
1575 /* Have to do this part _after_ outputting the logfile and the edr file */
1576 /* Gets written into the state at the beginning of next loop*/
1577 state->fep_state = lamnew;
1579 /* Print the remaining wall clock time for the run */
1580 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1582 if (shellfc)
1584 fprintf(stderr, "\n");
1586 print_time(stderr, walltime_accounting, step, ir, cr);
1589 /* Ion/water position swapping.
1590 * Not done in last step since trajectory writing happens before this call
1591 * in the MD loop and exchanges would be lost anyway. */
1592 bNeedRepartition = FALSE;
1593 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep && do_per_step(step, ir->swap->nstswap))
1595 bNeedRepartition =
1596 do_swapcoords(cr, step, t, ir, swap, wcycle, as_rvec_array(state->x.data()),
1597 state->box, MASTER(cr) && mdrunOptions.verbose, bRerunMD);
1599 if (bNeedRepartition && DOMAINDECOMP(cr))
1601 dd_collect_state(cr->dd, state, state_global);
1605 /* Replica exchange */
1606 bExchanged = FALSE;
1607 if (bDoReplEx)
1609 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1612 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1614 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1, state_global, *top_global, ir,
1615 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
1616 nrnb, wcycle, FALSE);
1617 shouldCheckNumberOfBondedInteractions = true;
1618 upd.setNumAtoms(state->natoms);
1621 bFirstStep = FALSE;
1622 bInitStep = FALSE;
1624 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1625 /* With all integrators, except VV, we need to retain the pressure
1626 * at the current step for coupling at the next step.
1628 if ((state->flags & (1U << estPRES_PREV))
1629 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1631 /* Store the pressure in t_state for pressure coupling
1632 * at the next MD step.
1634 copy_mat(pres, state->pres_prev);
1637 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1639 if ((membed != nullptr) && (!bLastStep))
1641 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1644 cycles = wallcycle_stop(wcycle, ewcSTEP);
1645 if (DOMAINDECOMP(cr) && wcycle)
1647 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1650 /* increase the MD step number */
1651 step++;
1652 step_rel++;
1654 resetHandler->resetCounters(step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb,
1655 fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1657 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1658 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1660 /* End of main MD loop */
1662 /* Closing TNG files can include compressing data. Therefore it is good to do that
1663 * before stopping the time measurements. */
1664 mdoutf_tng_close(outf);
1666 /* Stop measuring walltime */
1667 walltime_accounting_end_time(walltime_accounting);
1669 if (!thisRankHasDuty(cr, DUTY_PME))
1671 /* Tell the PME only node to finish */
1672 gmx_pme_send_finish(cr);
1675 if (MASTER(cr))
1677 if (ir->nstcalcenergy > 0)
1679 gmx::EnergyOutput::printAnnealingTemperatures(fplog, groups, &(ir->opts));
1680 energyOutput.printAverages(fplog, groups);
1683 done_mdoutf(outf);
1685 if (bPMETune)
1687 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1690 done_shellfc(fplog, shellfc, step_rel);
1692 if (useReplicaExchange && MASTER(cr))
1694 print_replica_exchange_statistics(fplog, repl_ex);
1697 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1699 global_stat_destroy(gstat);