2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019,2020, 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
52 #include "gromacs/awh/awh.h"
53 #include "gromacs/commandline/filenm.h"
54 #include "gromacs/domdec/collect.h"
55 #include "gromacs/domdec/dlbtiming.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/domdec/mdsetup.h"
60 #include "gromacs/domdec/partition.h"
61 #include "gromacs/essentialdynamics/edsam.h"
62 #include "gromacs/ewald/pme_load_balancing.h"
63 #include "gromacs/ewald/pme_pp.h"
64 #include "gromacs/fileio/trxio.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gmxlib/nrnb.h"
67 #include "gromacs/gpu_utils/gpu_utils.h"
68 #include "gromacs/listed_forces/manage_threading.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/utilities.h"
71 #include "gromacs/math/vec.h"
72 #include "gromacs/math/vectypes.h"
73 #include "gromacs/mdlib/checkpointhandler.h"
74 #include "gromacs/mdlib/compute_io.h"
75 #include "gromacs/mdlib/constr.h"
76 #include "gromacs/mdlib/ebin.h"
77 #include "gromacs/mdlib/enerdata_utils.h"
78 #include "gromacs/mdlib/energyoutput.h"
79 #include "gromacs/mdlib/expanded.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/md_support.h"
84 #include "gromacs/mdlib/mdatoms.h"
85 #include "gromacs/mdlib/mdoutf.h"
86 #include "gromacs/mdlib/membed.h"
87 #include "gromacs/mdlib/resethandler.h"
88 #include "gromacs/mdlib/sighandler.h"
89 #include "gromacs/mdlib/simulationsignal.h"
90 #include "gromacs/mdlib/stat.h"
91 #include "gromacs/mdlib/stophandler.h"
92 #include "gromacs/mdlib/tgroup.h"
93 #include "gromacs/mdlib/trajectory_writing.h"
94 #include "gromacs/mdlib/update.h"
95 #include "gromacs/mdlib/vcm.h"
96 #include "gromacs/mdlib/vsite.h"
97 #include "gromacs/mdrunutility/handlerestart.h"
98 #include "gromacs/mdrunutility/multisim.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 "legacysimulator.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,
155 const InteractionDefinitions& idef,
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.pbcType, 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, forceRec.pbcType, forceRec.bMolPBC, nullptr,
182 unshift_self(graph, globalState->box, globalState->x.rvec_array());
187 void gmx::LegacySimulator::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;
200 unsigned int force_flags;
201 tensor force_vir, shake_vir, total_vir, pres;
202 t_trxstatus* status = nullptr;
205 gmx_localtop_t top(top_global->ffparams);
206 PaddedHostVector<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);
230 "Note that it is planned that the command gmx mdrun -rerun will "
231 "be available in a different form in a future version of GROMACS, "
232 "e.g. gmx rerun -f.");
234 if (ir->efep != efepNO
235 && (mdAtoms->mdatoms()->nMassPerturbed > 0 || (constr && constr->havePerturbedConstraints())))
238 "Perturbed masses or constraints are not supported by rerun. "
239 "Either make a .tpr without mass and constraint perturbation, "
240 "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
244 gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
248 gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
252 gmx_fatal(FARGS, "AWH not supported by rerun.");
254 if (replExParams.exchangeInterval > 0)
256 gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
258 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
260 gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
264 gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
268 gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
270 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
271 [](int i) { return i != eannNO; }))
273 gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
276 /* Rerun can't work if an output file name is the same as the input file name.
277 * If this is the case, the user will get an error telling them what the issue is.
279 if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0
280 || strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
283 "When using mdrun -rerun, the name of the input trajectory file "
284 "%s cannot be identical to the name of an output file (whether "
285 "given explicitly with -o or -x, or by default)",
286 opt2fn("-rerun", nfile, fnm));
289 /* Settings for rerun */
291 ir->nstcalcenergy = 1;
292 int nstglobalcomm = 1;
293 const bool bNS = true;
295 ir->nstxout_compressed = 0;
296 SimulationGroups* groups = &top_global->groups;
297 if (ir->eI == eiMimic)
299 top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
302 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
303 const bool simulationsShareState = false;
304 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
305 mdModulesNotifier, ir, top_global, oenv, wcycle,
306 StartingBehavior::NewSimulation, simulationsShareState, ms);
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 stateInstance = std::make_unique<t_state>();
332 state = stateInstance.get();
333 dd_init_local_state(cr->dd, state_global, state);
335 /* Distribute the charge groups over the nodes from the master node */
336 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
337 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
338 nrnb, nullptr, FALSE);
339 shouldCheckNumberOfBondedInteractions = true;
343 state_change_natoms(state_global, state_global->natoms);
344 /* We need to allocate one element extra, since we might use
345 * (unaligned) 4-wide SIMD loads to access rvec entries.
347 f.resizeWithPadding(state_global->natoms);
348 /* Copy the pointer to the global state */
349 state = state_global;
351 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &graph, mdAtoms, constr, vsite, shellfc);
354 auto mdatoms = mdAtoms->mdatoms();
356 // NOTE: The global state is no longer used at this point.
357 // But state_global is still used as temporary storage space for writing
358 // the global state to file and potentially for replica exchange.
359 // (Global topology should persist.)
361 update_mdatoms(mdatoms, state->lambda[efptMASS]);
363 if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
365 doFreeEnergyPerturbation = true;
371 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
372 bool bSumEkinhOld = false;
373 t_vcm* vcm = nullptr;
374 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
375 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms, nrnb,
376 vcm, nullptr, enerd, force_vir, shake_vir, total_vir, pres, constr, &nullSignaller,
377 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
379 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
380 makeConstArrayRef(state->x), state->box,
381 &shouldCheckNumberOfBondedInteractions);
386 "starting md rerun '%s', reading coordinates from"
387 " input trajectory '%s'\n\n",
388 *(top_global->name), opt2fn("-rerun", nfile, fnm));
389 if (mdrunOptions.verbose)
392 "Calculated time to finish depends on nsteps from "
393 "run input file,\nwhich may not correspond to the time "
394 "needed to process input trajectory.\n\n");
396 fprintf(fplog, "\n");
399 walltime_accounting_start_time(walltime_accounting);
400 wallcycle_start(wcycle, ewcRUN);
401 print_start(fplog, cr, walltime_accounting, "mdrun");
403 /***********************************************************
407 ************************************************************/
413 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
419 isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
420 if (rerun_fr.natoms != top_global->natoms)
423 "Number of atoms in trajectory (%d) does not match the "
424 "run input file (%d)\n",
425 rerun_fr.natoms, top_global->natoms);
428 if (ir->pbcType != PbcType::No)
433 "Rerun trajectory frame step %" PRId64
435 "does not contain a box, while pbc is used",
436 rerun_fr.step, rerun_fr.time);
438 if (max_cutoff2(ir->pbcType, rerun_fr.box) < gmx::square(fr->rlist))
441 "Rerun trajectory frame step %" PRId64
443 "has too small box dimensions",
444 rerun_fr.step, rerun_fr.time);
452 "Rerun does not report kinetic energy, total energy, temperature, virial and "
457 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
460 if (ir->pbcType != PbcType::No)
462 /* Set the shift vectors.
463 * Necessary here when have a static box different from the tpr box.
465 calc_shifts(rerun_fr.box, fr->shift_vec);
468 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
469 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false, MASTER(cr),
470 ir->nstlist, mdrunOptions.reproducible, nstglobalcomm, mdrunOptions.maximumHoursToRun,
471 ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
473 // we don't do counter resetting in rerun - finish will always be valid
474 walltime_accounting_set_valid_finish(walltime_accounting);
476 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
478 step = ir->init_step;
481 /* and stop now if we should */
482 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
485 wallcycle_start(wcycle, ewcSTEP);
489 step = rerun_fr.step;
490 step_rel = step - ir->init_step;
501 if (ir->efep != efepNO && MASTER(cr))
503 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
508 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
509 if (constructVsites && DOMAINDECOMP(cr))
512 "Vsite recalculation with -rerun is not implemented with domain "
514 "use a single rank");
516 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top.idef, ir->delta_t,
520 isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
522 if (DOMAINDECOMP(cr))
524 /* Repartition the domain decomposition */
525 const bool bMasterState = true;
526 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
527 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
528 fr, vsite, constr, nrnb, wcycle, mdrunOptions.verbose);
529 shouldCheckNumberOfBondedInteractions = true;
534 energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
537 if (ir->efep != efepNO)
539 update_mdatoms(mdatoms, state->lambda[efptMASS]);
542 force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
543 | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
544 GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
548 /* Now is the time to relax the shells */
549 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
550 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
551 state->natoms, state->x.arrayRefWithPadding(),
552 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
553 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, graph,
554 shellfc, fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
558 /* The coordinates (x) are shifted (to get whole molecules)
560 * This is parallellized as well, and does communication too.
561 * Check comments in sim_util.c
564 gmx_edsam* ed = nullptr;
565 do_force(fplog, cr, ms, ir, awh, enforcedRotation, imdSession, pull_work, step, nrnb,
566 wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
567 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, graph,
568 fr, runScheduleWork, vsite, mu_tot, t, ed, GMX_FORCE_NS | force_flags,
569 ddBalanceRegionHandler);
572 /* Now we have the energies and forces corresponding to the
573 * coordinates at time t.
576 const bool isCheckpointingStep = false;
577 const bool doRerun = true;
578 const bool bSumEkinhOld = false;
579 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state,
580 state_global, observablesHistory, top_global, fr, outf,
581 energyOutput, ekind, f, isCheckpointingStep, doRerun,
582 isLastStep, mdrunOptions.writeConfout, bSumEkinhOld);
585 stopHandler->setSignal();
589 /* Need to unshift here */
590 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
593 if (vsite != nullptr)
595 wallcycle_start(wcycle, ewcVSITECONSTR);
596 if (graph != nullptr)
598 shift_self(graph, state->box, as_rvec_array(state->x.data()));
600 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t,
601 as_rvec_array(state->v.data()), top.idef.iparams, top.idef.il,
602 fr->pbcType, fr->bMolPBC, cr, state->box);
604 if (graph != nullptr)
606 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
608 wallcycle_stop(wcycle, ewcVSITECONSTR);
612 const bool doInterSimSignal = false;
613 const bool doIntraSimSignal = true;
614 bool bSumEkinhOld = false;
615 t_vcm* vcm = nullptr;
616 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
618 compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
619 makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
620 nrnb, vcm, wcycle, enerd, force_vir, shake_vir, total_vir, pres, constr,
621 &signaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
622 CGLO_GSTAT | CGLO_ENERGY
623 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
625 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
626 &top, makeConstArrayRef(state->x), state->box,
627 &shouldCheckNumberOfBondedInteractions);
630 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
631 the virial that should probably be addressed eventually. state->veta has better properies,
632 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
633 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
635 if (ir->efep != efepNO)
637 /* Sum up the foreign energy and dhdl terms for md and sd.
638 Currently done every step so that dhdl is correct in the .edr */
639 sum_dhdl(enerd, state->lambda, *ir->fepvals);
645 const bool bCalcEnerStep = true;
646 energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep, t,
647 mdatoms->tmass, enerd, state, ir->fepvals,
648 ir->expandedvals, state->box, shake_vir, force_vir,
649 total_vir, pres, ekind, mu_tot, constr);
651 const bool do_ene = true;
652 const bool do_log = true;
654 const bool do_dr = ir->nstdisreout != 0;
655 const bool do_or = ir->nstorireout != 0;
657 energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
658 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
659 do_log ? fplog : nullptr, step, t, fcd, awh);
661 if (do_per_step(step, ir->nstlog))
663 if (fflush(fplog) != 0)
665 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
670 /* Print the remaining wall clock time for the run */
671 if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
675 fprintf(stderr, "\n");
677 print_time(stderr, walltime_accounting, step, ir, cr);
680 /* Ion/water position swapping.
681 * Not done in last step since trajectory writing happens before this call
682 * in the MD loop and exchanges would be lost anyway. */
683 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep && do_per_step(step, ir->swap->nstswap))
685 const bool doRerun = true;
686 do_swapcoords(cr, step, t, ir, swap, wcycle, rerun_fr.x, rerun_fr.box,
687 MASTER(cr) && mdrunOptions.verbose, doRerun);
692 /* read next frame from input trajectory */
693 isLastStep = !read_next_frame(oenv, status, &rerun_fr);
698 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
701 cycles = wallcycle_stop(wcycle, ewcSTEP);
702 if (DOMAINDECOMP(cr) && wcycle)
704 dd_cycles_add(cr->dd, cycles, ddCyclStep);
709 /* increase the MD step number */
714 /* End of main MD loop */
716 /* Closing TNG files can include compressing data. Therefore it is good to do that
717 * before stopping the time measurements. */
718 mdoutf_tng_close(outf);
720 /* Stop measuring walltime */
721 walltime_accounting_end_time(walltime_accounting);
728 if (!thisRankHasDuty(cr, DUTY_PME))
730 /* Tell the PME only node to finish */
731 gmx_pme_send_finish(cr);
736 done_shellfc(fplog, shellfc, step_rel);
738 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);