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/compat/make_unique.h"
55 #include "gromacs/domdec/collect.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/partition.h"
60 #include "gromacs/essentialdynamics/edsam.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/ewald/pme-load-balancing.h"
63 #include "gromacs/fileio/trxio.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gmxlib/nrnb.h"
66 #include "gromacs/gpu_utils/gpu_utils.h"
67 #include "gromacs/imd/imd.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/expanded.h"
78 #include "gromacs/mdlib/force.h"
79 #include "gromacs/mdlib/force_flags.h"
80 #include "gromacs/mdlib/forcerec.h"
81 #include "gromacs/mdlib/md_support.h"
82 #include "gromacs/mdlib/mdatoms.h"
83 #include "gromacs/mdlib/mdebin.h"
84 #include "gromacs/mdlib/mdoutf.h"
85 #include "gromacs/mdlib/mdrun.h"
86 #include "gromacs/mdlib/mdsetup.h"
87 #include "gromacs/mdlib/membed.h"
88 #include "gromacs/mdlib/nb_verlet.h"
89 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
90 #include "gromacs/mdlib/ns.h"
91 #include "gromacs/mdlib/resethandler.h"
92 #include "gromacs/mdlib/shellfc.h"
93 #include "gromacs/mdlib/sighandler.h"
94 #include "gromacs/mdlib/sim_util.h"
95 #include "gromacs/mdlib/simulationsignal.h"
96 #include "gromacs/mdlib/stophandler.h"
97 #include "gromacs/mdlib/tgroup.h"
98 #include "gromacs/mdlib/trajectory_writing.h"
99 #include "gromacs/mdlib/update.h"
100 #include "gromacs/mdlib/vcm.h"
101 #include "gromacs/mdlib/vsite.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/observableshistory.h"
115 #include "gromacs/mdtypes/state.h"
116 #include "gromacs/mimic/MimicUtils.h"
117 #include "gromacs/pbcutil/mshift.h"
118 #include "gromacs/pbcutil/pbc.h"
119 #include "gromacs/pulling/pull.h"
120 #include "gromacs/swap/swapcoords.h"
121 #include "gromacs/timing/wallcycle.h"
122 #include "gromacs/timing/walltime_accounting.h"
123 #include "gromacs/topology/atoms.h"
124 #include "gromacs/topology/idef.h"
125 #include "gromacs/topology/mtop_util.h"
126 #include "gromacs/topology/topology.h"
127 #include "gromacs/trajectory/trajectoryframe.h"
128 #include "gromacs/utility/basedefinitions.h"
129 #include "gromacs/utility/cstringutil.h"
130 #include "gromacs/utility/fatalerror.h"
131 #include "gromacs/utility/logger.h"
132 #include "gromacs/utility/real.h"
133 #include "gromacs/utility/smalloc.h"
135 #include "integrator.h"
136 #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 gmx_mdoutf *outf = nullptr;
197 int64_t step, step_rel;
198 double t, lam0[efptNR];
199 bool isLastStep = false;
200 bool doFreeEnergyPerturbation = false;
202 tensor force_vir, shake_vir, total_vir, pres;
207 t_mdebin *mdebin = nullptr;
208 gmx_enerdata_t *enerd;
209 PaddedVector<gmx::RVec> f {};
210 gmx_global_stat_t gstat;
211 t_graph *graph = nullptr;
212 gmx_groups_t *groups;
213 gmx_ekindata_t *ekind;
214 gmx_shellfc_t *shellfc;
218 /* Domain decomposition could incorrectly miss a bonded
219 interaction, but checking for that requires a global
220 communication stage, which does not otherwise happen in DD
221 code. So we do that alongside the first global energy reduction
222 after a new DD is made. These variables handle whether the
223 check happens, and the result it returns. */
224 bool shouldCheckNumberOfBondedInteractions = false;
225 int totalNumberOfBondedInteractions = -1;
227 SimulationSignals signals;
228 // Most global communnication stages don't propagate mdrun
229 // signals, and will use this object to achieve that.
230 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
232 GMX_LOG(mdlog.info).asParagraph().
233 appendText("Note that it is planned that the command gmx mdrun -rerun will "
234 "be available in a different form in a future version of GROMACS, "
235 "e.g. gmx rerun -f.");
237 if (ir->efep != efepNO && (mdAtoms->mdatoms()->nMassPerturbed > 0 ||
238 (constr && constr->havePerturbedConstraints())))
240 gmx_fatal(FARGS, "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)
284 gmx_fatal(FARGS, "When using mdrun -rerun, the name of the input trajectory file "
285 "%s cannot be identical to the name of an output file (whether "
286 "given explicitly with -o or -x, or by default)",
287 opt2fn("-rerun", nfile, fnm));
290 /* Settings for rerun */
292 ir->nstcalcenergy = 1;
293 int nstglobalcomm = 1;
294 const bool bNS = true;
296 ir->nstxout_compressed = 0;
297 groups = &top_global->groups;
298 if (ir->eI == eiMimic)
300 top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
304 init_rerun(fplog, cr, outputProvider, ir, oenv, mdrunOptions,
305 state_global, lam0, nrnb, top_global,
306 nfile, fnm, &outf, &mdebin, wcycle);
308 /* Energy terms and groups */
310 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
313 /* Kinetic energy data */
315 init_ekindata(fplog, top_global, &(ir->opts), ekind);
316 /* Copy the cos acceleration to the groups struct */
317 ekind->cosacc.cos_accel = ir->cos_accel;
319 gstat = global_stat_init(ir);
321 /* Check for polarizable models and flexible constraints */
322 shellfc = init_shell_flexcon(fplog,
323 top_global, constr ? constr->numFlexibleConstraints() : 0,
324 ir->nstcalcenergy, DOMAINDECOMP(cr));
327 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
328 if ((io > 2000) && MASTER(cr))
331 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
336 // Local state only becomes valid now.
337 std::unique_ptr<t_state> stateInstance;
340 if (DOMAINDECOMP(cr))
342 top = dd_init_local_top(top_global);
344 stateInstance = compat::make_unique<t_state>();
345 state = stateInstance.get();
346 dd_init_local_state(cr->dd, state_global, state);
348 /* Distribute the charge groups over the nodes from the master node */
349 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
350 state_global, top_global, ir,
351 state, &f, mdAtoms, top, fr,
353 nrnb, nullptr, FALSE);
354 shouldCheckNumberOfBondedInteractions = true;
358 state_change_natoms(state_global, state_global->natoms);
359 /* We need to allocate one element extra, since we might use
360 * (unaligned) 4-wide SIMD loads to access rvec entries.
362 f.resizeWithPadding(state_global->natoms);
363 /* Copy the pointer to the global state */
364 state = state_global;
367 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
368 &graph, mdAtoms, constr, vsite, shellfc);
371 auto mdatoms = mdAtoms->mdatoms();
373 // NOTE: The global state is no longer used at this point.
374 // But state_global is still used as temporary storage space for writing
375 // the global state to file and potentially for replica exchange.
376 // (Global topology should persist.)
378 update_mdatoms(mdatoms, state->lambda[efptMASS]);
380 if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
382 doFreeEnergyPerturbation = true;
386 int cglo_flags = (CGLO_INITIALIZATION | CGLO_GSTAT |
387 (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
388 bool bSumEkinhOld = false;
389 t_vcm *vcm = nullptr;
390 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
391 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
392 constr, &nullSignaller, state->box,
393 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
395 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
396 top_global, top, state,
397 &shouldCheckNumberOfBondedInteractions);
401 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
402 " input trajectory '%s'\n\n",
403 *(top_global->name), opt2fn("-rerun", nfile, fnm));
404 if (mdrunOptions.verbose)
406 fprintf(stderr, "Calculated time to finish depends on nsteps from "
407 "run input file,\nwhich may not correspond to the time "
408 "needed to process input trajectory.\n\n");
410 fprintf(fplog, "\n");
413 walltime_accounting_start_time(walltime_accounting);
414 wallcycle_start(wcycle, ewcRUN);
415 print_start(fplog, cr, walltime_accounting, "mdrun");
417 /***********************************************************
421 ************************************************************/
425 GMX_LOG(mdlog.info).asParagraph().
426 appendText("Simulations has constraints. Rerun does not recalculate constraints.");
432 isLastStep = !read_first_frame(oenv, &status,
433 opt2fn("-rerun", nfile, fnm),
434 &rerun_fr, TRX_NEED_X);
435 if (rerun_fr.natoms != top_global->natoms)
438 "Number of atoms in trajectory (%d) does not match the "
439 "run input file (%d)\n",
440 rerun_fr.natoms, top_global->natoms);
443 if (ir->ePBC != epbcNONE)
447 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f "
448 "does not contain a box, while pbc is used",
449 rerun_fr.step, rerun_fr.time);
451 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
453 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f "
454 "has too small box dimensions", rerun_fr.step, rerun_fr.time);
459 GMX_LOG(mdlog.info).asParagraph().
460 appendText("Rerun does not report kinetic energy, total energy, temperature, virial and pressure.");
464 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
467 if (ir->ePBC != epbcNONE)
469 /* Set the shift vectors.
470 * Necessary here when have a static box different from the tpr box.
472 calc_shifts(rerun_fr.box, fr->shift_vec);
475 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
476 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false,
477 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstglobalcomm,
478 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
480 // we don't do counter resetting in rerun - finish will always be valid
481 walltime_accounting_set_valid_finish(walltime_accounting);
483 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion = (DOMAINDECOMP(cr) ? DdOpenBalanceRegionBeforeForceComputation::yes : DdOpenBalanceRegionBeforeForceComputation::no);
484 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion = (DOMAINDECOMP(cr) ? DdCloseBalanceRegionAfterForceComputation::yes : DdCloseBalanceRegionAfterForceComputation::no);
486 step = ir->init_step;
489 /* and stop now if we should */
490 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
493 wallcycle_start(wcycle, ewcSTEP);
497 step = rerun_fr.step;
498 step_rel = step - ir->init_step;
509 if (ir->efep != efepNO && MASTER(cr))
511 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
516 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
517 if (constructVsites && DOMAINDECOMP(cr))
519 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, "
520 "use a single rank");
522 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top->idef, ir->delta_t, *fr, graph);
525 isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
527 if (DOMAINDECOMP(cr))
529 /* Repartition the domain decomposition */
530 const bool bMasterState = true;
531 dd_partition_system(fplog, mdlog, step, cr,
532 bMasterState, nstglobalcomm,
533 state_global, top_global, ir,
534 state, &f, mdAtoms, top, fr,
537 mdrunOptions.verbose);
538 shouldCheckNumberOfBondedInteractions = true;
543 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
546 if (ir->efep != efepNO)
548 update_mdatoms(mdatoms, state->lambda[efptMASS]);
551 force_flags = (GMX_FORCE_STATECHANGED |
552 GMX_FORCE_DYNAMICBOX |
553 GMX_FORCE_ALLFORCES |
554 GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
556 (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
560 /* Now is the time to relax the shells */
561 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
562 enforcedRotation, step,
563 ir, bNS, force_flags, top,
565 state, f.arrayRefWithPadding(), force_vir, mdatoms,
566 nrnb, wcycle, graph, groups,
567 shellfc, fr, ppForceWorkload, t, mu_tot,
569 ddOpenBalanceRegion, ddCloseBalanceRegion);
573 /* The coordinates (x) are shifted (to get whole molecules)
575 * This is parallellized as well, and does communication too.
576 * Check comments in sim_util.c
579 gmx_edsam *ed = nullptr;
580 do_force(fplog, cr, ms, ir, awh, enforcedRotation,
581 step, nrnb, wcycle, top, groups,
582 state->box, state->x.arrayRefWithPadding(), &state->hist,
583 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
584 state->lambda, graph,
585 fr, ppForceWorkload, vsite, mu_tot, t, ed,
586 GMX_FORCE_NS | force_flags,
587 ddOpenBalanceRegion, ddCloseBalanceRegion);
590 /* Now we have the energies and forces corresponding to the
591 * coordinates at time t.
594 const bool isCheckpointingStep = false;
595 const bool doRerun = true;
596 const bool bSumEkinhOld = false;
597 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
598 ir, state, state_global, observablesHistory,
600 outf, mdebin, ekind, f,
601 isCheckpointingStep, doRerun, isLastStep,
602 mdrunOptions.writeConfout,
606 stopHandler->setSignal();
610 /* Need to unshift here */
611 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
614 if (vsite != nullptr)
616 wallcycle_start(wcycle, ewcVSITECONSTR);
617 if (graph != nullptr)
619 shift_self(graph, state->box, as_rvec_array(state->x.data()));
621 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
622 top->idef.iparams, top->idef.il,
623 fr->ePBC, fr->bMolPBC, cr, state->box);
625 if (graph != nullptr)
627 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
629 wallcycle_stop(wcycle, ewcVSITECONSTR);
633 const bool doInterSimSignal = false;
634 const bool doIntraSimSignal = true;
635 bool bSumEkinhOld = false;
636 t_vcm *vcm = nullptr;
637 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
639 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
640 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
643 &totalNumberOfBondedInteractions, &bSumEkinhOld,
644 CGLO_GSTAT | CGLO_ENERGY
645 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
647 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
648 top_global, top, state,
649 &shouldCheckNumberOfBondedInteractions);
652 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
653 the virial that should probably be addressed eventually. state->veta has better properies,
654 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
655 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
657 if (ir->efep != efepNO)
659 /* Sum up the foreign energy and dhdl terms for md and sd.
660 Currently done every step so that dhdl is correct in the .edr */
661 sum_dhdl(enerd, state->lambda, ir->fepvals);
667 const bool bCalcEnerStep = true;
668 upd_mdebin(mdebin, doFreeEnergyPerturbation, bCalcEnerStep,
669 t, mdatoms->tmass, enerd, state,
670 ir->fepvals, ir->expandedvals, rerun_fr.box,
671 shake_vir, force_vir, total_vir, pres,
672 ekind, mu_tot, constr);
674 const bool do_ene = true;
675 const bool do_log = true;
677 const bool do_dr = ir->nstdisreout != 0;
678 const bool do_or = ir->nstorireout != 0;
680 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
682 eprNORMAL, mdebin, fcd, groups, &(ir->opts), awh);
684 if (do_per_step(step, ir->nstlog))
686 if (fflush(fplog) != 0)
688 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
693 /* Print the remaining wall clock time for the run */
694 if (isMasterSimMasterRank(ms, cr) &&
695 (mdrunOptions.verbose || gmx_got_usr_signal()))
699 fprintf(stderr, "\n");
701 print_time(stderr, walltime_accounting, step, ir, cr);
704 /* Ion/water position swapping.
705 * Not done in last step since trajectory writing happens before this call
706 * in the MD loop and exchanges would be lost anyway. */
707 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep &&
708 do_per_step(step, ir->swap->nstswap))
710 const bool doRerun = true;
711 do_swapcoords(cr, step, t, ir, wcycle,
714 MASTER(cr) && mdrunOptions.verbose,
720 /* read next frame from input trajectory */
721 isLastStep = !read_next_frame(oenv, status, &rerun_fr);
726 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
729 cycles = wallcycle_stop(wcycle, ewcSTEP);
730 if (DOMAINDECOMP(cr) && wcycle)
732 dd_cycles_add(cr->dd, cycles, ddCyclStep);
737 /* increase the MD step number */
742 /* End of main MD loop */
744 /* Closing TNG files can include compressing data. Therefore it is good to do that
745 * before stopping the time measurements. */
746 mdoutf_tng_close(outf);
748 /* Stop measuring walltime */
749 walltime_accounting_end_time(walltime_accounting);
756 if (!thisRankHasDuty(cr, DUTY_PME))
758 /* Tell the PME only node to finish */
759 gmx_pme_send_finish(cr);
765 done_shellfc(fplog, shellfc, step_rel);
767 // Clean up swapcoords
768 if (ir->eSwapCoords != eswapNO)
770 finish_swapcoords(ir->swap);
773 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
775 destroy_enerdata(enerd);