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