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/handlerestart.h"
100 #include "gromacs/mdrunutility/multisim.h"
101 #include "gromacs/mdrunutility/printtime.h"
102 #include "gromacs/mdtypes/awh_history.h"
103 #include "gromacs/mdtypes/awh_params.h"
104 #include "gromacs/mdtypes/commrec.h"
105 #include "gromacs/mdtypes/df_history.h"
106 #include "gromacs/mdtypes/energyhistory.h"
107 #include "gromacs/mdtypes/fcdata.h"
108 #include "gromacs/mdtypes/forcerec.h"
109 #include "gromacs/mdtypes/group.h"
110 #include "gromacs/mdtypes/inputrec.h"
111 #include "gromacs/mdtypes/interaction_const.h"
112 #include "gromacs/mdtypes/md_enums.h"
113 #include "gromacs/mdtypes/mdatom.h"
114 #include "gromacs/mdtypes/mdrunoptions.h"
115 #include "gromacs/mdtypes/observableshistory.h"
116 #include "gromacs/mdtypes/state.h"
117 #include "gromacs/mimic/utilities.h"
118 #include "gromacs/pbcutil/mshift.h"
119 #include "gromacs/pbcutil/pbc.h"
120 #include "gromacs/pulling/pull.h"
121 #include "gromacs/swap/swapcoords.h"
122 #include "gromacs/timing/wallcycle.h"
123 #include "gromacs/timing/walltime_accounting.h"
124 #include "gromacs/topology/atoms.h"
125 #include "gromacs/topology/idef.h"
126 #include "gromacs/topology/mtop_util.h"
127 #include "gromacs/topology/topology.h"
128 #include "gromacs/trajectory/trajectoryframe.h"
129 #include "gromacs/utility/basedefinitions.h"
130 #include "gromacs/utility/cstringutil.h"
131 #include "gromacs/utility/fatalerror.h"
132 #include "gromacs/utility/logger.h"
133 #include "gromacs/utility/real.h"
134 #include "gromacs/utility/smalloc.h"
136 #include "legacysimulator.h"
137 #include "replicaexchange.h"
140 using gmx::SimulationSignaller;
142 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
144 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
145 * \param[in,out] globalState The global state container
146 * \param[in] constructVsites When true, vsite coordinates are constructed
147 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
148 * \param[in] idef Topology parameters, used for constructing vsites
149 * \param[in] timeStep Time step, used for constructing vsites
150 * \param[in] forceRec Force record, used for constructing vsites
151 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
153 static void prepareRerunState(const t_trxframe& rerunFrame,
154 t_state* globalState,
155 bool constructVsites,
156 const gmx_vsite_t* vsite,
159 const t_forcerec& forceRec,
162 auto x = makeArrayRef(globalState->x);
163 auto rerunX = arrayRefFromArray(reinterpret_cast<gmx::RVec*>(rerunFrame.x), globalState->natoms);
164 std::copy(rerunX.begin(), rerunX.end(), x.begin());
165 copy_mat(rerunFrame.box, globalState->box);
169 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
173 /* Following is necessary because the graph may get out of sync
174 * with the coordinates if we only have every N'th coordinate set
176 mk_mshift(nullptr, graph, forceRec.ePBC, globalState->box, globalState->x.rvec_array());
177 shift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
179 construct_vsites(vsite, globalState->x.rvec_array(), timeStep, globalState->v.rvec_array(),
180 idef.iparams, idef.il, forceRec.ePBC, forceRec.bMolPBC, nullptr,
184 unshift_self(graph, globalState->box, globalState->x.rvec_array());
189 void gmx::LegacySimulator::do_rerun()
191 // TODO Historically, the EM and MD "integrators" used different
192 // names for the t_inputrec *parameter, but these must have the
193 // same name, now that it's a member of a struct. We use this ir
194 // alias to avoid a large ripple of nearly useless changes.
195 // t_inputrec is being replaced by IMdpOptionsProvider, so this
196 // will go away eventually.
197 t_inputrec* ir = inputrec;
198 int64_t step, step_rel;
199 double t, lam0[efptNR];
200 bool isLastStep = false;
201 bool doFreeEnergyPerturbation = false;
202 unsigned int force_flags;
203 tensor force_vir, shake_vir, total_vir, pres;
204 t_trxstatus* status = nullptr;
208 PaddedHostVector<gmx::RVec> f{};
209 gmx_global_stat_t gstat;
210 t_graph* graph = nullptr;
211 gmx_shellfc_t* shellfc;
215 /* Domain decomposition could incorrectly miss a bonded
216 interaction, but checking for that requires a global
217 communication stage, which does not otherwise happen in DD
218 code. So we do that alongside the first global energy reduction
219 after a new DD is made. These variables handle whether the
220 check happens, and the result it returns. */
221 bool shouldCheckNumberOfBondedInteractions = false;
222 int totalNumberOfBondedInteractions = -1;
224 SimulationSignals signals;
225 // Most global communnication stages don't propagate mdrun
226 // signals, and will use this object to achieve that.
227 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
232 "Note that it is planned that the command gmx mdrun -rerun will "
233 "be available in a different form in a future version of GROMACS, "
234 "e.g. gmx rerun -f.");
236 if (ir->efep != efepNO
237 && (mdAtoms->mdatoms()->nMassPerturbed > 0 || (constr && constr->havePerturbedConstraints())))
240 "Perturbed masses or constraints are not supported by rerun. "
241 "Either make a .tpr without mass and constraint perturbation, "
242 "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
246 gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
250 gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
254 gmx_fatal(FARGS, "AWH not supported by rerun.");
256 if (replExParams.exchangeInterval > 0)
258 gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
260 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
262 gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
266 gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
270 gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
272 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
273 [](int i) { return i != eannNO; }))
275 gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
278 /* Rerun can't work if an output file name is the same as the input file name.
279 * If this is the case, the user will get an error telling them what the issue is.
281 if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0
282 || strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
285 "When using mdrun -rerun, the name of the input trajectory file "
286 "%s cannot be identical to the name of an output file (whether "
287 "given explicitly with -o or -x, or by default)",
288 opt2fn("-rerun", nfile, fnm));
291 /* Settings for rerun */
293 ir->nstcalcenergy = 1;
294 int nstglobalcomm = 1;
295 const bool bNS = true;
297 ir->nstxout_compressed = 0;
298 SimulationGroups* groups = &top_global->groups;
299 if (ir->eI == eiMimic)
301 top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
304 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
305 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier,
306 ir, top_global, oenv, wcycle, StartingBehavior::NewSimulation);
307 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
308 mdoutf_get_fp_dhdl(outf), true, StartingBehavior::NewSimulation,
311 gstat = global_stat_init(ir);
313 /* Check for polarizable models and flexible constraints */
314 shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
315 ir->nstcalcenergy, DOMAINDECOMP(cr));
318 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
319 if ((io > 2000) && MASTER(cr))
321 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
325 // Local state only becomes valid now.
326 std::unique_ptr<t_state> stateInstance;
329 if (DOMAINDECOMP(cr))
331 dd_init_local_top(*top_global, &top);
333 stateInstance = std::make_unique<t_state>();
334 state = stateInstance.get();
335 dd_init_local_state(cr->dd, state_global, state);
337 /* Distribute the charge groups over the nodes from the master node */
338 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
339 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
340 nrnb, nullptr, FALSE);
341 shouldCheckNumberOfBondedInteractions = true;
345 state_change_natoms(state_global, state_global->natoms);
346 /* We need to allocate one element extra, since we might use
347 * (unaligned) 4-wide SIMD loads to access rvec entries.
349 f.resizeWithPadding(state_global->natoms);
350 /* Copy the pointer to the global state */
351 state = state_global;
353 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &graph, mdAtoms, constr, vsite, shellfc);
356 auto mdatoms = mdAtoms->mdatoms();
358 // NOTE: The global state is no longer used at this point.
359 // But state_global is still used as temporary storage space for writing
360 // the global state to file and potentially for replica exchange.
361 // (Global topology should persist.)
363 update_mdatoms(mdatoms, state->lambda[efptMASS]);
365 if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
367 doFreeEnergyPerturbation = true;
373 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
374 bool bSumEkinhOld = false;
375 t_vcm* vcm = nullptr;
376 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
377 state->box, state->lambda[efptVDW], mdatoms, nrnb, vcm, nullptr, enerd,
378 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
379 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
381 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
382 state->x.rvec_array(), state->box,
383 &shouldCheckNumberOfBondedInteractions);
388 "starting md rerun '%s', reading coordinates from"
389 " input trajectory '%s'\n\n",
390 *(top_global->name), opt2fn("-rerun", nfile, fnm));
391 if (mdrunOptions.verbose)
394 "Calculated time to finish depends on nsteps from "
395 "run input file,\nwhich may not correspond to the time "
396 "needed to process input trajectory.\n\n");
398 fprintf(fplog, "\n");
401 walltime_accounting_start_time(walltime_accounting);
402 wallcycle_start(wcycle, ewcRUN);
403 print_start(fplog, cr, walltime_accounting, "mdrun");
405 /***********************************************************
409 ************************************************************/
415 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
421 isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
422 if (rerun_fr.natoms != top_global->natoms)
425 "Number of atoms in trajectory (%d) does not match the "
426 "run input file (%d)\n",
427 rerun_fr.natoms, top_global->natoms);
430 if (ir->ePBC != epbcNONE)
435 "Rerun trajectory frame step %" PRId64
437 "does not contain a box, while pbc is used",
438 rerun_fr.step, rerun_fr.time);
440 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
443 "Rerun trajectory frame step %" PRId64
445 "has too small box dimensions",
446 rerun_fr.step, rerun_fr.time);
454 "Rerun does not report kinetic energy, total energy, temperature, virial and "
459 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
462 if (ir->ePBC != epbcNONE)
464 /* Set the shift vectors.
465 * Necessary here when have a static box different from the tpr box.
467 calc_shifts(rerun_fr.box, fr->shift_vec);
470 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
471 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false, MASTER(cr),
472 ir->nstlist, mdrunOptions.reproducible, nstglobalcomm, mdrunOptions.maximumHoursToRun,
473 ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
475 // we don't do counter resetting in rerun - finish will always be valid
476 walltime_accounting_set_valid_finish(walltime_accounting);
478 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
480 step = ir->init_step;
483 /* and stop now if we should */
484 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
487 wallcycle_start(wcycle, ewcSTEP);
491 step = rerun_fr.step;
492 step_rel = step - ir->init_step;
503 if (ir->efep != efepNO && MASTER(cr))
505 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
510 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
511 if (constructVsites && DOMAINDECOMP(cr))
514 "Vsite recalculation with -rerun is not implemented with domain "
516 "use a single rank");
518 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top.idef, ir->delta_t,
522 isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
524 if (DOMAINDECOMP(cr))
526 /* Repartition the domain decomposition */
527 const bool bMasterState = true;
528 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
529 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
530 fr, vsite, constr, nrnb, wcycle, mdrunOptions.verbose);
531 shouldCheckNumberOfBondedInteractions = true;
536 energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
539 if (ir->efep != efepNO)
541 update_mdatoms(mdatoms, state->lambda[efptMASS]);
544 force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
545 | (GMX_GPU ? GMX_FORCE_VIRIAL : 0) | // TODO: Get rid of this once #2649 is solved
546 GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
550 /* Now is the time to relax the shells */
551 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
552 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
553 state->natoms, state->x.arrayRefWithPadding(),
554 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
555 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, graph,
556 shellfc, fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
560 /* The coordinates (x) are shifted (to get whole molecules)
562 * This is parallellized as well, and does communication too.
563 * Check comments in sim_util.c
566 gmx_edsam* ed = nullptr;
567 do_force(fplog, cr, ms, ir, awh, enforcedRotation, imdSession, pull_work, step, nrnb,
568 wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
569 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, graph,
570 fr, runScheduleWork, vsite, mu_tot, t, ed, GMX_FORCE_NS | force_flags,
571 ddBalanceRegionHandler);
574 /* Now we have the energies and forces corresponding to the
575 * coordinates at time t.
578 const bool isCheckpointingStep = false;
579 const bool doRerun = true;
580 const bool bSumEkinhOld = false;
581 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state,
582 state_global, observablesHistory, top_global, fr, outf,
583 energyOutput, ekind, f, isCheckpointingStep, doRerun,
584 isLastStep, mdrunOptions.writeConfout, bSumEkinhOld);
587 stopHandler->setSignal();
591 /* Need to unshift here */
592 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
595 if (vsite != nullptr)
597 wallcycle_start(wcycle, ewcVSITECONSTR);
598 if (graph != nullptr)
600 shift_self(graph, state->box, as_rvec_array(state->x.data()));
602 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t,
603 as_rvec_array(state->v.data()), top.idef.iparams, top.idef.il,
604 fr->ePBC, fr->bMolPBC, cr, state->box);
606 if (graph != nullptr)
608 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
610 wallcycle_stop(wcycle, ewcVSITECONSTR);
614 const bool doInterSimSignal = false;
615 const bool doIntraSimSignal = true;
616 bool bSumEkinhOld = false;
617 t_vcm* vcm = nullptr;
618 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
620 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
621 state->box, state->lambda[efptVDW], mdatoms, nrnb, vcm, wcycle, enerd,
622 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &signaller,
623 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
624 CGLO_GSTAT | CGLO_ENERGY
625 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
627 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
628 &top, state->x.rvec_array(), state->box,
629 &shouldCheckNumberOfBondedInteractions);
632 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
633 the virial that should probably be addressed eventually. state->veta has better properies,
634 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
635 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
637 if (ir->efep != efepNO)
639 /* Sum up the foreign energy and dhdl terms for md and sd.
640 Currently done every step so that dhdl is correct in the .edr */
641 sum_dhdl(enerd, state->lambda, *ir->fepvals);
647 const bool bCalcEnerStep = true;
648 energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep, t,
649 mdatoms->tmass, enerd, state, ir->fepvals,
650 ir->expandedvals, state->box, shake_vir, force_vir,
651 total_vir, pres, ekind, mu_tot, constr);
653 const bool do_ene = true;
654 const bool do_log = true;
656 const bool do_dr = ir->nstdisreout != 0;
657 const bool do_or = ir->nstorireout != 0;
659 energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
660 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
661 do_log ? fplog : nullptr, step, t, fcd, awh);
663 if (do_per_step(step, ir->nstlog))
665 if (fflush(fplog) != 0)
667 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
672 /* Print the remaining wall clock time for the run */
673 if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
677 fprintf(stderr, "\n");
679 print_time(stderr, walltime_accounting, step, ir, cr);
682 /* Ion/water position swapping.
683 * Not done in last step since trajectory writing happens before this call
684 * in the MD loop and exchanges would be lost anyway. */
685 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep && do_per_step(step, ir->swap->nstswap))
687 const bool doRerun = true;
688 do_swapcoords(cr, step, t, ir, swap, wcycle, rerun_fr.x, rerun_fr.box,
689 MASTER(cr) && mdrunOptions.verbose, doRerun);
694 /* read next frame from input trajectory */
695 isLastStep = !read_next_frame(oenv, status, &rerun_fr);
700 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
703 cycles = wallcycle_stop(wcycle, ewcSTEP);
704 if (DOMAINDECOMP(cr) && wcycle)
706 dd_cycles_add(cr->dd, cycles, ddCyclStep);
711 /* increase the MD step number */
716 /* End of main MD loop */
718 /* Closing TNG files can include compressing data. Therefore it is good to do that
719 * before stopping the time measurements. */
720 mdoutf_tng_close(outf);
722 /* Stop measuring walltime */
723 walltime_accounting_end_time(walltime_accounting);
730 if (!thisRankHasDuty(cr, DUTY_PME))
732 /* Tell the PME only node to finish */
733 gmx_pme_send_finish(cr);
738 done_shellfc(fplog, shellfc, step_rel);
740 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);