2c59d9635dd0aabff1bd7bc851ecc2d737b194ad
[alexxy/gromacs.git] / src / gromacs / mdrun / rerun.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief Implements the loop for simulation reruns
38  *
39  * \author Pascal Merz <pascal.merz@colorado.edu>
40  * \ingroup module_mdrun
41  */
42 #include "gmxpre.h"
43
44 #include <cinttypes>
45 #include <cmath>
46 #include <cstdio>
47 #include <cstdlib>
48
49 #include <algorithm>
50 #include <memory>
51
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/listed_forces.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/forcerec.h"
106 #include "gromacs/mdtypes/group.h"
107 #include "gromacs/mdtypes/inputrec.h"
108 #include "gromacs/mdtypes/interaction_const.h"
109 #include "gromacs/mdtypes/md_enums.h"
110 #include "gromacs/mdtypes/mdatom.h"
111 #include "gromacs/mdtypes/mdrunoptions.h"
112 #include "gromacs/mdtypes/observableshistory.h"
113 #include "gromacs/mdtypes/state.h"
114 #include "gromacs/mimic/utilities.h"
115 #include "gromacs/pbcutil/pbc.h"
116 #include "gromacs/pulling/pull.h"
117 #include "gromacs/swap/swapcoords.h"
118 #include "gromacs/timing/wallcycle.h"
119 #include "gromacs/timing/walltime_accounting.h"
120 #include "gromacs/topology/atoms.h"
121 #include "gromacs/topology/idef.h"
122 #include "gromacs/topology/mtop_util.h"
123 #include "gromacs/topology/topology.h"
124 #include "gromacs/trajectory/trajectoryframe.h"
125 #include "gromacs/utility/basedefinitions.h"
126 #include "gromacs/utility/cstringutil.h"
127 #include "gromacs/utility/fatalerror.h"
128 #include "gromacs/utility/logger.h"
129 #include "gromacs/utility/real.h"
130 #include "gromacs/utility/smalloc.h"
131
132 #include "legacysimulator.h"
133 #include "replicaexchange.h"
134 #include "shellfc.h"
135
136 using gmx::SimulationSignaller;
137 using gmx::VirtualSitesHandler;
138
139 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
140  *
141  * \param[in]     rerunFrame      The trajectory frame to compute energy/forces for
142  * \param[in,out] globalState     The global state container
143  * \param[in]     constructVsites When true, vsite coordinates are constructed
144  * \param[in]     vsite           Vsite setup, can be nullptr when \p constructVsites = false
145  * \param[in]     timeStep        Time step, used for constructing vsites
146  */
147 static void prepareRerunState(const t_trxframe&          rerunFrame,
148                               t_state*                   globalState,
149                               bool                       constructVsites,
150                               const VirtualSitesHandler* vsite,
151                               double                     timeStep)
152 {
153     auto x      = makeArrayRef(globalState->x);
154     auto rerunX = arrayRefFromArray(reinterpret_cast<gmx::RVec*>(rerunFrame.x), globalState->natoms);
155     std::copy(rerunX.begin(), rerunX.end(), x.begin());
156     copy_mat(rerunFrame.box, globalState->box);
157
158     if (constructVsites)
159     {
160         GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
161
162         vsite->construct(globalState->x, timeStep, globalState->v, globalState->box);
163     }
164 }
165
166 void gmx::LegacySimulator::do_rerun()
167 {
168     // TODO Historically, the EM and MD "integrators" used different
169     // names for the t_inputrec *parameter, but these must have the
170     // same name, now that it's a member of a struct. We use this ir
171     // alias to avoid a large ripple of nearly useless changes.
172     // t_inputrec is being replaced by IMdpOptionsProvider, so this
173     // will go away eventually.
174     t_inputrec*                 ir = inputrec;
175     int64_t                     step, step_rel;
176     double                      t, lam0[efptNR];
177     bool                        isLastStep               = false;
178     bool                        doFreeEnergyPerturbation = false;
179     unsigned int                force_flags;
180     tensor                      force_vir, shake_vir, total_vir, pres;
181     t_trxstatus*                status = nullptr;
182     rvec                        mu_tot;
183     t_trxframe                  rerun_fr;
184     gmx_localtop_t              top(top_global->ffparams);
185     PaddedHostVector<gmx::RVec> f{};
186     gmx_global_stat_t           gstat;
187     gmx_shellfc_t*              shellfc;
188
189     double cycles;
190
191     /* Domain decomposition could incorrectly miss a bonded
192        interaction, but checking for that requires a global
193        communication stage, which does not otherwise happen in DD
194        code. So we do that alongside the first global energy reduction
195        after a new DD is made. These variables handle whether the
196        check happens, and the result it returns. */
197     bool shouldCheckNumberOfBondedInteractions = false;
198     int  totalNumberOfBondedInteractions       = -1;
199
200     SimulationSignals signals;
201     // Most global communnication stages don't propagate mdrun
202     // signals, and will use this object to achieve that.
203     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
204
205     GMX_LOG(mdlog.info)
206             .asParagraph()
207             .appendText(
208                     "Note that it is planned that the command gmx mdrun -rerun will "
209                     "be available in a different form in a future version of GROMACS, "
210                     "e.g. gmx rerun -f.");
211
212     if (ir->efep != efepNO
213         && (mdAtoms->mdatoms()->nMassPerturbed > 0 || (constr && constr->havePerturbedConstraints())))
214     {
215         gmx_fatal(FARGS,
216                   "Perturbed masses or constraints are not supported by rerun. "
217                   "Either make a .tpr without mass and constraint perturbation, "
218                   "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
219     }
220     if (ir->bExpanded)
221     {
222         gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
223     }
224     if (ir->bSimTemp)
225     {
226         gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
227     }
228     if (ir->bDoAwh)
229     {
230         gmx_fatal(FARGS, "AWH not supported by rerun.");
231     }
232     if (replExParams.exchangeInterval > 0)
233     {
234         gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
235     }
236     if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
237     {
238         gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
239     }
240     if (ir->bIMD)
241     {
242         gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
243     }
244     if (isMultiSim(ms))
245     {
246         gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
247     }
248     if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
249                     [](int i) { return i != eannNO; }))
250     {
251         gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
252     }
253
254     /* Rerun can't work if an output file name is the same as the input file name.
255      * If this is the case, the user will get an error telling them what the issue is.
256      */
257     if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0
258         || strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
259     {
260         gmx_fatal(FARGS,
261                   "When using mdrun -rerun, the name of the input trajectory file "
262                   "%s cannot be identical to the name of an output file (whether "
263                   "given explicitly with -o or -x, or by default)",
264                   opt2fn("-rerun", nfile, fnm));
265     }
266
267     /* Settings for rerun */
268     ir->nstlist              = 1;
269     ir->nstcalcenergy        = 1;
270     int        nstglobalcomm = 1;
271     const bool bNS           = true;
272
273     ir->nstxout_compressed         = 0;
274     const SimulationGroups* groups = &top_global->groups;
275     if (ir->eI == eiMimic)
276     {
277         auto nonConstGlobalTopology                          = const_cast<gmx_mtop_t*>(top_global);
278         nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(*top_global);
279     }
280
281     initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
282     const bool        simulationsShareState = false;
283     gmx_mdoutf*       outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
284                                    mdModulesNotifier, ir, top_global, oenv, wcycle,
285                                    StartingBehavior::NewSimulation, simulationsShareState, ms);
286     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
287                                    mdoutf_get_fp_dhdl(outf), true, StartingBehavior::NewSimulation,
288                                    mdModulesNotifier);
289
290     gstat = global_stat_init(ir);
291
292     /* Check for polarizable models and flexible constraints */
293     shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
294                                  ir->nstcalcenergy, DOMAINDECOMP(cr));
295
296     {
297         double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
298         if ((io > 2000) && MASTER(cr))
299         {
300             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
301         }
302     }
303
304     // Local state only becomes valid now.
305     std::unique_ptr<t_state> stateInstance;
306     t_state*                 state;
307
308     if (DOMAINDECOMP(cr))
309     {
310         stateInstance = std::make_unique<t_state>();
311         state         = stateInstance.get();
312         dd_init_local_state(cr->dd, state_global, state);
313
314         /* Distribute the charge groups over the nodes from the master node */
315         dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
316                             imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
317                             nrnb, nullptr, FALSE);
318         shouldCheckNumberOfBondedInteractions = true;
319     }
320     else
321     {
322         state_change_natoms(state_global, state_global->natoms);
323         /* Copy the pointer to the global state */
324         state = state_global;
325
326         mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
327     }
328
329     auto mdatoms = mdAtoms->mdatoms();
330
331     // NOTE: The global state is no longer used at this point.
332     // But state_global is still used as temporary storage space for writing
333     // the global state to file and potentially for replica exchange.
334     // (Global topology should persist.)
335
336     update_mdatoms(mdatoms, state->lambda[efptMASS]);
337
338     if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
339     {
340         doFreeEnergyPerturbation = true;
341     }
342
343     {
344         int cglo_flags =
345                 (CGLO_GSTAT
346                  | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
347         bool   bSumEkinhOld = false;
348         t_vcm* vcm          = nullptr;
349         compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
350                         makeConstArrayRef(state->v), state->box, mdatoms, nrnb, vcm, nullptr, enerd,
351                         force_vir, shake_vir, total_vir, pres, constr, &nullSignaller, state->box,
352                         &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
353     }
354     checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
355                                     makeConstArrayRef(state->x), state->box,
356                                     &shouldCheckNumberOfBondedInteractions);
357
358     if (MASTER(cr))
359     {
360         fprintf(stderr,
361                 "starting md rerun '%s', reading coordinates from"
362                 " input trajectory '%s'\n\n",
363                 *(top_global->name), opt2fn("-rerun", nfile, fnm));
364         if (mdrunOptions.verbose)
365         {
366             fprintf(stderr,
367                     "Calculated time to finish depends on nsteps from "
368                     "run input file,\nwhich may not correspond to the time "
369                     "needed to process input trajectory.\n\n");
370         }
371         fprintf(fplog, "\n");
372     }
373
374     walltime_accounting_start_time(walltime_accounting);
375     wallcycle_start(wcycle, ewcRUN);
376     print_start(fplog, cr, walltime_accounting, "mdrun");
377
378     /***********************************************************
379      *
380      *             Loop over MD steps
381      *
382      ************************************************************/
383
384     if (constr)
385     {
386         GMX_LOG(mdlog.info)
387                 .asParagraph()
388                 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
389     }
390
391     rerun_fr.natoms = 0;
392     if (MASTER(cr))
393     {
394         isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
395         if (rerun_fr.natoms != top_global->natoms)
396         {
397             gmx_fatal(FARGS,
398                       "Number of atoms in trajectory (%d) does not match the "
399                       "run input file (%d)\n",
400                       rerun_fr.natoms, top_global->natoms);
401         }
402
403         if (ir->pbcType != PbcType::No)
404         {
405             if (!rerun_fr.bBox)
406             {
407                 gmx_fatal(FARGS,
408                           "Rerun trajectory frame step %" PRId64
409                           " time %f "
410                           "does not contain a box, while pbc is used",
411                           rerun_fr.step, rerun_fr.time);
412             }
413             if (max_cutoff2(ir->pbcType, rerun_fr.box) < gmx::square(fr->rlist))
414             {
415                 gmx_fatal(FARGS,
416                           "Rerun trajectory frame step %" PRId64
417                           " time %f "
418                           "has too small box dimensions",
419                           rerun_fr.step, rerun_fr.time);
420             }
421         }
422     }
423
424     GMX_LOG(mdlog.info)
425             .asParagraph()
426             .appendText(
427                     "Rerun does not report kinetic energy, total energy, temperature, virial and "
428                     "pressure.");
429
430     if (PAR(cr))
431     {
432         rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
433     }
434
435     if (ir->pbcType != PbcType::No)
436     {
437         /* Set the shift vectors.
438          * Necessary here when have a static box different from the tpr box.
439          */
440         calc_shifts(rerun_fr.box, fr->shift_vec);
441     }
442
443     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
444             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false, MASTER(cr),
445             ir->nstlist, mdrunOptions.reproducible, nstglobalcomm, mdrunOptions.maximumHoursToRun,
446             ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
447
448     // we don't do counter resetting in rerun - finish will always be valid
449     walltime_accounting_set_valid_finish(walltime_accounting);
450
451     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
452
453     step     = ir->init_step;
454     step_rel = 0;
455
456     /* and stop now if we should */
457     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
458     while (!isLastStep)
459     {
460         wallcycle_start(wcycle, ewcSTEP);
461
462         if (rerun_fr.bStep)
463         {
464             step     = rerun_fr.step;
465             step_rel = step - ir->init_step;
466         }
467         if (rerun_fr.bTime)
468         {
469             t = rerun_fr.time;
470         }
471         else
472         {
473             t = step;
474         }
475
476         if (ir->efep != efepNO && MASTER(cr))
477         {
478             setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
479         }
480
481         if (MASTER(cr))
482         {
483             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
484             if (constructVsites && DOMAINDECOMP(cr))
485             {
486                 gmx_fatal(FARGS,
487                           "Vsite recalculation with -rerun is not implemented with domain "
488                           "decomposition, "
489                           "use a single rank");
490             }
491             prepareRerunState(rerun_fr, state_global, constructVsites, vsite, ir->delta_t);
492         }
493
494         isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
495
496         if (DOMAINDECOMP(cr))
497         {
498             /* Repartition the domain decomposition */
499             const bool bMasterState = true;
500             dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
501                                 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
502                                 fr, vsite, constr, nrnb, wcycle, mdrunOptions.verbose);
503             shouldCheckNumberOfBondedInteractions = true;
504         }
505
506         if (MASTER(cr))
507         {
508             EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
509         }
510
511         if (ir->efep != efepNO)
512         {
513             update_mdatoms(mdatoms, state->lambda[efptMASS]);
514         }
515
516         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
517                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
518                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
519
520         if (shellfc)
521         {
522             /* Now is the time to relax the shells */
523             relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
524                                 imdSession, pull_work, bNS, force_flags, &top, constr, enerd,
525                                 state->natoms, state->x.arrayRefWithPadding(),
526                                 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
527                                 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, shellfc,
528                                 fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
529         }
530         else
531         {
532             /* The coordinates (x) are shifted (to get whole molecules)
533              * in do_force.
534              * This is parallellized as well, and does communication too.
535              * Check comments in sim_util.c
536              */
537             Awh*       awh = nullptr;
538             gmx_edsam* ed  = nullptr;
539             do_force(fplog, cr, ms, ir, awh, enforcedRotation, imdSession, pull_work, step, nrnb,
540                      wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
541                      f.arrayRefWithPadding(), force_vir, mdatoms, enerd, state->lambda, fr, runScheduleWork,
542                      vsite, mu_tot, t, ed, GMX_FORCE_NS | force_flags, ddBalanceRegionHandler);
543         }
544
545         /* Now we have the energies and forces corresponding to the
546          * coordinates at time t.
547          */
548         {
549             const bool isCheckpointingStep = false;
550             const bool doRerun             = true;
551             const bool bSumEkinhOld        = false;
552             do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state,
553                                      state_global, observablesHistory, top_global, fr, outf,
554                                      energyOutput, ekind, f, isCheckpointingStep, doRerun,
555                                      isLastStep, mdrunOptions.writeConfout, bSumEkinhOld);
556         }
557
558         stopHandler->setSignal();
559
560         if (vsite != nullptr)
561         {
562             wallcycle_start(wcycle, ewcVSITECONSTR);
563             vsite->construct(state->x, ir->delta_t, state->v, state->box);
564             wallcycle_stop(wcycle, ewcVSITECONSTR);
565         }
566
567         {
568             const bool          doInterSimSignal = false;
569             const bool          doIntraSimSignal = true;
570             bool                bSumEkinhOld     = false;
571             t_vcm*              vcm              = nullptr;
572             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
573
574             compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
575                             makeConstArrayRef(state->v), state->box, mdatoms, nrnb, vcm, wcycle,
576                             enerd, force_vir, shake_vir, total_vir, pres, constr, &signaller,
577                             state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
578                             CGLO_GSTAT | CGLO_ENERGY
579                                     | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
580                                                                              : 0));
581             checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
582                                             &top, makeConstArrayRef(state->x), state->box,
583                                             &shouldCheckNumberOfBondedInteractions);
584         }
585
586         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
587            the virial that should probably be addressed eventually. state->veta has better properies,
588            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
589            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
590
591         if (ir->efep != efepNO)
592         {
593             /* Sum up the foreign energy and dhdl terms for md and sd.
594                Currently done every step so that dhdl is correct in the .edr */
595             sum_dhdl(enerd, state->lambda, *ir->fepvals);
596         }
597
598         /* Output stuff */
599         if (MASTER(cr))
600         {
601             const bool bCalcEnerStep = true;
602             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep, t,
603                                              mdatoms->tmass, enerd, state, ir->fepvals,
604                                              ir->expandedvals, state->box, shake_vir, force_vir,
605                                              total_vir, pres, ekind, mu_tot, constr);
606
607             const bool do_ene = true;
608             const bool do_log = true;
609             Awh*       awh    = nullptr;
610             const bool do_dr  = ir->nstdisreout != 0;
611             const bool do_or  = ir->nstorireout != 0;
612
613             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
614             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
615                                                do_log ? fplog : nullptr, step, t,
616                                                &fr->listedForces->fcdata(), awh);
617
618             if (do_per_step(step, ir->nstlog))
619             {
620                 if (fflush(fplog) != 0)
621                 {
622                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
623                 }
624             }
625         }
626
627         /* Print the remaining wall clock time for the run */
628         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
629         {
630             if (shellfc)
631             {
632                 fprintf(stderr, "\n");
633             }
634             print_time(stderr, walltime_accounting, step, ir, cr);
635         }
636
637         /* Ion/water position swapping.
638          * Not done in last step since trajectory writing happens before this call
639          * in the MD loop and exchanges would be lost anyway. */
640         if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep && do_per_step(step, ir->swap->nstswap))
641         {
642             const bool doRerun = true;
643             do_swapcoords(cr, step, t, ir, swap, wcycle, rerun_fr.x, rerun_fr.box,
644                           MASTER(cr) && mdrunOptions.verbose, doRerun);
645         }
646
647         if (MASTER(cr))
648         {
649             /* read next frame from input trajectory */
650             isLastStep = !read_next_frame(oenv, status, &rerun_fr);
651         }
652
653         if (PAR(cr))
654         {
655             rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
656         }
657
658         cycles = wallcycle_stop(wcycle, ewcSTEP);
659         if (DOMAINDECOMP(cr) && wcycle)
660         {
661             dd_cycles_add(cr->dd, cycles, ddCyclStep);
662         }
663
664         if (!rerun_fr.bStep)
665         {
666             /* increase the MD step number */
667             step++;
668             step_rel++;
669         }
670     }
671     /* End of main MD loop */
672
673     /* Closing TNG files can include compressing data. Therefore it is good to do that
674      * before stopping the time measurements. */
675     mdoutf_tng_close(outf);
676
677     /* Stop measuring walltime */
678     walltime_accounting_end_time(walltime_accounting);
679
680     if (MASTER(cr))
681     {
682         close_trx(status);
683     }
684
685     if (!thisRankHasDuty(cr, DUTY_PME))
686     {
687         /* Tell the PME only node to finish */
688         gmx_pme_send_finish(cr);
689     }
690
691     done_mdoutf(outf);
692
693     done_shellfc(fplog, shellfc, step_rel);
694
695     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
696 }