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