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