2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * \brief Implements the loop for simulation reruns
39 * \author Pascal Merz <pascal.merz@colorado.edu>
40 * \ingroup module_mdrun
54 #include "gromacs/awh/awh.h"
55 #include "gromacs/commandline/filenm.h"
56 #include "gromacs/domdec/collect.h"
57 #include "gromacs/domdec/dlbtiming.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_network.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/domdec/mdsetup.h"
62 #include "gromacs/domdec/partition.h"
63 #include "gromacs/essentialdynamics/edsam.h"
64 #include "gromacs/ewald/pme.h"
65 #include "gromacs/ewald/pme_load_balancing.h"
66 #include "gromacs/fileio/trxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/gpu_utils/gpu_utils.h"
70 #include "gromacs/listed_forces/manage_threading.h"
71 #include "gromacs/math/functions.h"
72 #include "gromacs/math/utilities.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/math/vectypes.h"
75 #include "gromacs/mdlib/checkpointhandler.h"
76 #include "gromacs/mdlib/compute_io.h"
77 #include "gromacs/mdlib/constr.h"
78 #include "gromacs/mdlib/ebin.h"
79 #include "gromacs/mdlib/enerdata_utils.h"
80 #include "gromacs/mdlib/energyoutput.h"
81 #include "gromacs/mdlib/expanded.h"
82 #include "gromacs/mdlib/force.h"
83 #include "gromacs/mdlib/force_flags.h"
84 #include "gromacs/mdlib/forcerec.h"
85 #include "gromacs/mdlib/md_support.h"
86 #include "gromacs/mdlib/mdatoms.h"
87 #include "gromacs/mdlib/mdoutf.h"
88 #include "gromacs/mdlib/membed.h"
89 #include "gromacs/mdlib/resethandler.h"
90 #include "gromacs/mdlib/sighandler.h"
91 #include "gromacs/mdlib/simulationsignal.h"
92 #include "gromacs/mdlib/stat.h"
93 #include "gromacs/mdlib/stophandler.h"
94 #include "gromacs/mdlib/tgroup.h"
95 #include "gromacs/mdlib/trajectory_writing.h"
96 #include "gromacs/mdlib/update.h"
97 #include "gromacs/mdlib/vcm.h"
98 #include "gromacs/mdlib/vsite.h"
99 #include "gromacs/mdrunutility/printtime.h"
100 #include "gromacs/mdtypes/awh_history.h"
101 #include "gromacs/mdtypes/awh_params.h"
102 #include "gromacs/mdtypes/commrec.h"
103 #include "gromacs/mdtypes/df_history.h"
104 #include "gromacs/mdtypes/energyhistory.h"
105 #include "gromacs/mdtypes/fcdata.h"
106 #include "gromacs/mdtypes/forcerec.h"
107 #include "gromacs/mdtypes/group.h"
108 #include "gromacs/mdtypes/inputrec.h"
109 #include "gromacs/mdtypes/interaction_const.h"
110 #include "gromacs/mdtypes/md_enums.h"
111 #include "gromacs/mdtypes/mdatom.h"
112 #include "gromacs/mdtypes/mdrunoptions.h"
113 #include "gromacs/mdtypes/observableshistory.h"
114 #include "gromacs/mdtypes/state.h"
115 #include "gromacs/mimic/utilities.h"
116 #include "gromacs/pbcutil/mshift.h"
117 #include "gromacs/pbcutil/pbc.h"
118 #include "gromacs/pulling/pull.h"
119 #include "gromacs/swap/swapcoords.h"
120 #include "gromacs/timing/wallcycle.h"
121 #include "gromacs/timing/walltime_accounting.h"
122 #include "gromacs/topology/atoms.h"
123 #include "gromacs/topology/idef.h"
124 #include "gromacs/topology/mtop_util.h"
125 #include "gromacs/topology/topology.h"
126 #include "gromacs/trajectory/trajectoryframe.h"
127 #include "gromacs/utility/basedefinitions.h"
128 #include "gromacs/utility/cstringutil.h"
129 #include "gromacs/utility/fatalerror.h"
130 #include "gromacs/utility/logger.h"
131 #include "gromacs/utility/real.h"
132 #include "gromacs/utility/smalloc.h"
134 #include "integrator.h"
135 #include "replicaexchange.h"
138 using gmx::SimulationSignaller;
140 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
142 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
143 * \param[in,out] globalState The global state container
144 * \param[in] constructVsites When true, vsite coordinates are constructed
145 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
146 * \param[in] idef Topology parameters, used for constructing vsites
147 * \param[in] timeStep Time step, used for constructing vsites
148 * \param[in] forceRec Force record, used for constructing vsites
149 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
151 static void prepareRerunState(const t_trxframe &rerunFrame,
152 t_state *globalState,
153 bool constructVsites,
154 const gmx_vsite_t *vsite,
157 const t_forcerec &forceRec,
160 auto x = makeArrayRef(globalState->x);
161 auto rerunX = arrayRefFromArray(reinterpret_cast<gmx::RVec *>(rerunFrame.x), globalState->natoms);
162 std::copy(rerunX.begin(), rerunX.end(), x.begin());
163 copy_mat(rerunFrame.box, globalState->box);
167 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
171 /* Following is necessary because the graph may get out of sync
172 * with the coordinates if we only have every N'th coordinate set
174 mk_mshift(nullptr, graph, forceRec.ePBC, globalState->box, globalState->x.rvec_array());
175 shift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
177 construct_vsites(vsite, globalState->x.rvec_array(), timeStep, globalState->v.rvec_array(),
178 idef.iparams, idef.il,
179 forceRec.ePBC, forceRec.bMolPBC, nullptr, globalState->box);
182 unshift_self(graph, globalState->box, globalState->x.rvec_array());
187 void gmx::Integrator::do_rerun()
189 // TODO Historically, the EM and MD "integrators" used different
190 // names for the t_inputrec *parameter, but these must have the
191 // same name, now that it's a member of a struct. We use this ir
192 // alias to avoid a large ripple of nearly useless changes.
193 // t_inputrec is being replaced by IMdpOptionsProvider, so this
194 // will go away eventually.
195 t_inputrec *ir = inputrec;
196 int64_t step, step_rel;
197 double t, lam0[efptNR];
198 bool isLastStep = false;
199 bool doFreeEnergyPerturbation = false;
201 tensor force_vir, shake_vir, total_vir, pres;
206 PaddedVector<gmx::RVec> f {};
207 gmx_global_stat_t gstat;
208 t_graph *graph = nullptr;
209 gmx_shellfc_t *shellfc;
213 /* Domain decomposition could incorrectly miss a bonded
214 interaction, but checking for that requires a global
215 communication stage, which does not otherwise happen in DD
216 code. So we do that alongside the first global energy reduction
217 after a new DD is made. These variables handle whether the
218 check happens, and the result it returns. */
219 bool shouldCheckNumberOfBondedInteractions = false;
220 int totalNumberOfBondedInteractions = -1;
222 SimulationSignals signals;
223 // Most global communnication stages don't propagate mdrun
224 // signals, and will use this object to achieve that.
225 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
227 GMX_LOG(mdlog.info).asParagraph().
228 appendText("Note that it is planned that the command gmx mdrun -rerun will "
229 "be available in a different form in a future version of GROMACS, "
230 "e.g. gmx rerun -f.");
232 if (ir->efep != efepNO && (mdAtoms->mdatoms()->nMassPerturbed > 0 ||
233 (constr && constr->havePerturbedConstraints())))
235 gmx_fatal(FARGS, "Perturbed masses or constraints are not supported by rerun. "
236 "Either make a .tpr without mass and constraint perturbation, "
237 "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
241 gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
245 gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
249 gmx_fatal(FARGS, "AWH not supported by rerun.");
251 if (replExParams.exchangeInterval > 0)
253 gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
255 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
257 gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
261 gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
265 gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
267 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
268 [](int i){return i != eannNO; }))
270 gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
273 /* Rerun can't work if an output file name is the same as the input file name.
274 * If this is the case, the user will get an error telling them what the issue is.
276 if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0 ||
277 strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
279 gmx_fatal(FARGS, "When using mdrun -rerun, the name of the input trajectory file "
280 "%s cannot be identical to the name of an output file (whether "
281 "given explicitly with -o or -x, or by default)",
282 opt2fn("-rerun", nfile, fnm));
285 /* Settings for rerun */
287 ir->nstcalcenergy = 1;
288 int nstglobalcomm = 1;
289 const bool bNS = true;
291 ir->nstxout_compressed = 0;
292 SimulationGroups *groups = &top_global->groups;
293 if (ir->eI == eiMimic)
295 top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
298 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
300 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, oenv, wcycle);
301 gmx::EnergyOutput energyOutput;
302 energyOutput.prepare(mdoutf_get_fp_ene(outf), top_global, ir, pull_work, mdoutf_get_fp_dhdl(outf), true);
304 /* Kinetic energy data */
305 std::unique_ptr<gmx_ekindata_t> eKinData = std::make_unique<gmx_ekindata_t>();
306 gmx_ekindata_t *ekind = eKinData.get();
307 init_ekindata(fplog, top_global, &(ir->opts), ekind);
308 /* Copy the cos acceleration to the groups struct */
309 ekind->cosacc.cos_accel = ir->cos_accel;
311 gstat = global_stat_init(ir);
313 /* Check for polarizable models and flexible constraints */
314 shellfc = init_shell_flexcon(fplog,
315 top_global, constr ? constr->numFlexibleConstraints() : 0,
316 ir->nstcalcenergy, DOMAINDECOMP(cr));
319 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
320 if ((io > 2000) && MASTER(cr))
323 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
328 // Local state only becomes valid now.
329 std::unique_ptr<t_state> stateInstance;
332 if (DOMAINDECOMP(cr))
334 dd_init_local_top(*top_global, &top);
336 stateInstance = std::make_unique<t_state>();
337 state = stateInstance.get();
338 dd_init_local_state(cr->dd, state_global, state);
340 /* Distribute the charge groups over the nodes from the master node */
341 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
342 state_global, *top_global, ir, imdSession,
344 state, &f, mdAtoms, &top, fr,
346 nrnb, nullptr, FALSE);
347 shouldCheckNumberOfBondedInteractions = true;
351 state_change_natoms(state_global, state_global->natoms);
352 /* We need to allocate one element extra, since we might use
353 * (unaligned) 4-wide SIMD loads to access rvec entries.
355 f.resizeWithPadding(state_global->natoms);
356 /* Copy the pointer to the global state */
357 state = state_global;
359 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr,
360 &graph, mdAtoms, constr, vsite, shellfc);
363 auto mdatoms = mdAtoms->mdatoms();
365 // NOTE: The global state is no longer used at this point.
366 // But state_global is still used as temporary storage space for writing
367 // the global state to file and potentially for replica exchange.
368 // (Global topology should persist.)
370 update_mdatoms(mdatoms, state->lambda[efptMASS]);
372 if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
374 doFreeEnergyPerturbation = true;
378 int cglo_flags = (CGLO_INITIALIZATION | CGLO_GSTAT |
379 (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
380 bool bSumEkinhOld = false;
381 t_vcm *vcm = nullptr;
382 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
383 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
384 constr, &nullSignaller, state->box,
385 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
387 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
388 top_global, &top, state,
389 &shouldCheckNumberOfBondedInteractions);
393 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
394 " input trajectory '%s'\n\n",
395 *(top_global->name), opt2fn("-rerun", nfile, fnm));
396 if (mdrunOptions.verbose)
398 fprintf(stderr, "Calculated time to finish depends on nsteps from "
399 "run input file,\nwhich may not correspond to the time "
400 "needed to process input trajectory.\n\n");
402 fprintf(fplog, "\n");
405 walltime_accounting_start_time(walltime_accounting);
406 wallcycle_start(wcycle, ewcRUN);
407 print_start(fplog, cr, walltime_accounting, "mdrun");
409 /***********************************************************
413 ************************************************************/
417 GMX_LOG(mdlog.info).asParagraph().
418 appendText("Simulations has constraints. Rerun does not recalculate constraints.");
424 isLastStep = !read_first_frame(oenv, &status,
425 opt2fn("-rerun", nfile, fnm),
426 &rerun_fr, TRX_NEED_X);
427 if (rerun_fr.natoms != top_global->natoms)
430 "Number of atoms in trajectory (%d) does not match the "
431 "run input file (%d)\n",
432 rerun_fr.natoms, top_global->natoms);
435 if (ir->ePBC != epbcNONE)
439 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f "
440 "does not contain a box, while pbc is used",
441 rerun_fr.step, rerun_fr.time);
443 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
445 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f "
446 "has too small box dimensions", rerun_fr.step, rerun_fr.time);
451 GMX_LOG(mdlog.info).asParagraph().
452 appendText("Rerun does not report kinetic energy, total energy, temperature, virial and pressure.");
456 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
459 if (ir->ePBC != epbcNONE)
461 /* Set the shift vectors.
462 * Necessary here when have a static box different from the tpr box.
464 calc_shifts(rerun_fr.box, fr->shift_vec);
467 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
468 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false,
469 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstglobalcomm,
470 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
472 // we don't do counter resetting in rerun - finish will always be valid
473 walltime_accounting_set_valid_finish(walltime_accounting);
475 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
477 step = ir->init_step;
480 /* and stop now if we should */
481 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
484 wallcycle_start(wcycle, ewcSTEP);
488 step = rerun_fr.step;
489 step_rel = step - ir->init_step;
500 if (ir->efep != efepNO && MASTER(cr))
502 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
507 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
508 if (constructVsites && DOMAINDECOMP(cr))
510 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, "
511 "use a single rank");
513 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top.idef, ir->delta_t, *fr, graph);
516 isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
518 if (DOMAINDECOMP(cr))
520 /* Repartition the domain decomposition */
521 const bool bMasterState = true;
522 dd_partition_system(fplog, mdlog, step, cr,
523 bMasterState, nstglobalcomm,
524 state_global, *top_global, ir, imdSession,
526 state, &f, mdAtoms, &top, fr,
529 mdrunOptions.verbose);
530 shouldCheckNumberOfBondedInteractions = true;
535 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
538 if (ir->efep != efepNO)
540 update_mdatoms(mdatoms, state->lambda[efptMASS]);
543 force_flags = (GMX_FORCE_STATECHANGED |
544 GMX_FORCE_DYNAMICBOX |
545 GMX_FORCE_ALLFORCES |
546 (GMX_GPU ? GMX_FORCE_VIRIAL : 0) | // TODO: Get rid of this once #2649 is solved
548 (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
552 /* Now is the time to relax the shells */
553 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
554 enforcedRotation, step,
555 ir, imdSession, pull_work, bNS, force_flags, &top,
557 state, f.arrayRefWithPadding(), force_vir, mdatoms,
559 shellfc, fr, ppForceWorkload, t, mu_tot,
561 ddBalanceRegionHandler);
565 /* The coordinates (x) are shifted (to get whole molecules)
567 * This is parallellized as well, and does communication too.
568 * Check comments in sim_util.c
571 gmx_edsam *ed = nullptr;
572 do_force(fplog, cr, ms, ir, awh, enforcedRotation, imdSession,
574 step, nrnb, wcycle, &top,
575 state->box, state->x.arrayRefWithPadding(), &state->hist,
576 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
577 state->lambda, graph,
578 fr, ppForceWorkload, vsite, mu_tot, t, ed,
579 GMX_FORCE_NS | force_flags,
580 ddBalanceRegionHandler);
583 /* Now we have the energies and forces corresponding to the
584 * coordinates at time t.
587 const bool isCheckpointingStep = false;
588 const bool doRerun = true;
589 const bool bSumEkinhOld = false;
590 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
591 ir, state, state_global, observablesHistory,
593 outf, energyOutput, ekind, f,
594 isCheckpointingStep, doRerun, isLastStep,
595 mdrunOptions.writeConfout,
599 stopHandler->setSignal();
603 /* Need to unshift here */
604 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
607 if (vsite != nullptr)
609 wallcycle_start(wcycle, ewcVSITECONSTR);
610 if (graph != nullptr)
612 shift_self(graph, state->box, as_rvec_array(state->x.data()));
614 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
615 top.idef.iparams, top.idef.il,
616 fr->ePBC, fr->bMolPBC, cr, state->box);
618 if (graph != nullptr)
620 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
622 wallcycle_stop(wcycle, ewcVSITECONSTR);
626 const bool doInterSimSignal = false;
627 const bool doIntraSimSignal = true;
628 bool bSumEkinhOld = false;
629 t_vcm *vcm = nullptr;
630 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
632 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
633 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
636 &totalNumberOfBondedInteractions, &bSumEkinhOld,
637 CGLO_GSTAT | CGLO_ENERGY
638 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
640 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
641 top_global, &top, state,
642 &shouldCheckNumberOfBondedInteractions);
645 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
646 the virial that should probably be addressed eventually. state->veta has better properies,
647 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
648 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
650 if (ir->efep != efepNO)
652 /* Sum up the foreign energy and dhdl terms for md and sd.
653 Currently done every step so that dhdl is correct in the .edr */
654 sum_dhdl(enerd, state->lambda, ir->fepvals);
660 const bool bCalcEnerStep = true;
661 energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep,
662 t, mdatoms->tmass, enerd, state,
663 ir->fepvals, ir->expandedvals, state->box,
664 shake_vir, force_vir, total_vir, pres,
665 ekind, mu_tot, constr);
667 const bool do_ene = true;
668 const bool do_log = true;
670 const bool do_dr = ir->nstdisreout != 0;
671 const bool do_or = ir->nstorireout != 0;
673 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
674 do_log ? fplog : nullptr,
676 eprNORMAL, fcd, groups, &(ir->opts), awh);
678 if (do_per_step(step, ir->nstlog))
680 if (fflush(fplog) != 0)
682 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
687 /* Print the remaining wall clock time for the run */
688 if (isMasterSimMasterRank(ms, cr) &&
689 (mdrunOptions.verbose || gmx_got_usr_signal()))
693 fprintf(stderr, "\n");
695 print_time(stderr, walltime_accounting, step, ir, cr);
698 /* Ion/water position swapping.
699 * Not done in last step since trajectory writing happens before this call
700 * in the MD loop and exchanges would be lost anyway. */
701 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep &&
702 do_per_step(step, ir->swap->nstswap))
704 const bool doRerun = true;
705 do_swapcoords(cr, step, t, ir, swap, wcycle,
708 MASTER(cr) && mdrunOptions.verbose,
714 /* read next frame from input trajectory */
715 isLastStep = !read_next_frame(oenv, status, &rerun_fr);
720 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
723 cycles = wallcycle_stop(wcycle, ewcSTEP);
724 if (DOMAINDECOMP(cr) && wcycle)
726 dd_cycles_add(cr->dd, cycles, ddCyclStep);
731 /* increase the MD step number */
736 /* End of main MD loop */
738 /* Closing TNG files can include compressing data. Therefore it is good to do that
739 * before stopping the time measurements. */
740 mdoutf_tng_close(outf);
742 /* Stop measuring walltime */
743 walltime_accounting_end_time(walltime_accounting);
750 if (!thisRankHasDuty(cr, DUTY_PME))
752 /* Tell the PME only node to finish */
753 gmx_pme_send_finish(cr);
758 done_shellfc(fplog, shellfc, step_rel);
760 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);