8230157ded92d002fd10ea995d2f99484102e7b0
[alexxy/gromacs.git] / src / gromacs / mdrun / rerun.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  *
37  * \brief Implements the loop for simulation reruns
38  *
39  * \author Pascal Merz <pascal.merz@colorado.edu>
40  * \ingroup module_mdrun
41  */
42 #include "gmxpre.h"
43
44 #include <cinttypes>
45 #include <cmath>
46 #include <cstdio>
47 #include <cstdlib>
48
49 #include <algorithm>
50 #include <memory>
51
52 #include "gromacs/awh/awh.h"
53 #include "gromacs/commandline/filenm.h"
54 #include "gromacs/domdec/collect.h"
55 #include "gromacs/domdec/dlbtiming.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/domdec/mdsetup.h"
60 #include "gromacs/domdec/partition.h"
61 #include "gromacs/essentialdynamics/edsam.h"
62 #include "gromacs/ewald/pme_load_balancing.h"
63 #include "gromacs/ewald/pme_pp.h"
64 #include "gromacs/fileio/trxio.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gmxlib/nrnb.h"
67 #include "gromacs/gpu_utils/gpu_utils.h"
68 #include "gromacs/listed_forces/listed_forces.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/utilities.h"
71 #include "gromacs/math/vec.h"
72 #include "gromacs/math/vectypes.h"
73 #include "gromacs/mdlib/checkpointhandler.h"
74 #include "gromacs/mdlib/compute_io.h"
75 #include "gromacs/mdlib/constr.h"
76 #include "gromacs/mdlib/ebin.h"
77 #include "gromacs/mdlib/enerdata_utils.h"
78 #include "gromacs/mdlib/energyoutput.h"
79 #include "gromacs/mdlib/expanded.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/md_support.h"
84 #include "gromacs/mdlib/mdatoms.h"
85 #include "gromacs/mdlib/mdoutf.h"
86 #include "gromacs/mdlib/membed.h"
87 #include "gromacs/mdlib/resethandler.h"
88 #include "gromacs/mdlib/sighandler.h"
89 #include "gromacs/mdlib/simulationsignal.h"
90 #include "gromacs/mdlib/stat.h"
91 #include "gromacs/mdlib/stophandler.h"
92 #include "gromacs/mdlib/tgroup.h"
93 #include "gromacs/mdlib/trajectory_writing.h"
94 #include "gromacs/mdlib/update.h"
95 #include "gromacs/mdlib/vcm.h"
96 #include "gromacs/mdlib/vsite.h"
97 #include "gromacs/mdrunutility/handlerestart.h"
98 #include "gromacs/mdrunutility/multisim.h"
99 #include "gromacs/mdrunutility/printtime.h"
100 #include "gromacs/mdtypes/awh_history.h"
101 #include "gromacs/mdtypes/awh_params.h"
102 #include "gromacs/mdtypes/commrec.h"
103 #include "gromacs/mdtypes/df_history.h"
104 #include "gromacs/mdtypes/energyhistory.h"
105 #include "gromacs/mdtypes/forcebuffers.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     ForceBuffers      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     int*                fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
282     gmx::ArrayRef<real> lambda    = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
283     initialize_lambdas(fplog, *ir, MASTER(cr), fep_state, lambda, lam0);
284     const bool        simulationsShareState = false;
285     gmx_mdoutf*       outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
286                                    mdModulesNotifier, ir, top_global, oenv, wcycle,
287                                    StartingBehavior::NewSimulation, simulationsShareState, ms);
288     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
289                                    mdoutf_get_fp_dhdl(outf), true, StartingBehavior::NewSimulation,
290                                    mdModulesNotifier);
291
292     gstat = global_stat_init(ir);
293
294     /* Check for polarizable models and flexible constraints */
295     shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
296                                  ir->nstcalcenergy, DOMAINDECOMP(cr));
297
298     {
299         double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
300         if ((io > 2000) && MASTER(cr))
301         {
302             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
303         }
304     }
305
306     // Local state only becomes valid now.
307     std::unique_ptr<t_state> stateInstance;
308     t_state*                 state;
309
310     if (DOMAINDECOMP(cr))
311     {
312         stateInstance = std::make_unique<t_state>();
313         state         = stateInstance.get();
314         dd_init_local_state(cr->dd, state_global, state);
315
316         /* Distribute the charge groups over the nodes from the master node */
317         dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
318                             imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
319                             nrnb, nullptr, FALSE);
320         shouldCheckNumberOfBondedInteractions = true;
321     }
322     else
323     {
324         state_change_natoms(state_global, state_global->natoms);
325         /* Copy the pointer to the global state */
326         state = state_global;
327
328         mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
329     }
330
331     auto mdatoms = mdAtoms->mdatoms();
332
333     // NOTE: The global state is no longer used at this point.
334     // But state_global is still used as temporary storage space for writing
335     // the global state to file and potentially for replica exchange.
336     // (Global topology should persist.)
337
338     update_mdatoms(mdatoms, state->lambda[efptMASS]);
339
340     if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
341     {
342         doFreeEnergyPerturbation = true;
343     }
344
345     {
346         int cglo_flags =
347                 (CGLO_GSTAT
348                  | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
349         bool   bSumEkinhOld = false;
350         t_vcm* vcm          = nullptr;
351         compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
352                         makeConstArrayRef(state->v), state->box, mdatoms, nrnb, vcm, nullptr, enerd,
353                         force_vir, shake_vir, total_vir, pres, constr, &nullSignaller, state->box,
354                         &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
355     }
356     checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
357                                     makeConstArrayRef(state->x), state->box,
358                                     &shouldCheckNumberOfBondedInteractions);
359
360     if (MASTER(cr))
361     {
362         fprintf(stderr,
363                 "starting md rerun '%s', reading coordinates from"
364                 " input trajectory '%s'\n\n",
365                 *(top_global->name), opt2fn("-rerun", nfile, fnm));
366         if (mdrunOptions.verbose)
367         {
368             fprintf(stderr,
369                     "Calculated time to finish depends on nsteps from "
370                     "run input file,\nwhich may not correspond to the time "
371                     "needed to process input trajectory.\n\n");
372         }
373         fprintf(fplog, "\n");
374     }
375
376     walltime_accounting_start_time(walltime_accounting);
377     wallcycle_start(wcycle, ewcRUN);
378     print_start(fplog, cr, walltime_accounting, "mdrun");
379
380     /***********************************************************
381      *
382      *             Loop over MD steps
383      *
384      ************************************************************/
385
386     if (constr)
387     {
388         GMX_LOG(mdlog.info)
389                 .asParagraph()
390                 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
391     }
392
393     rerun_fr.natoms = 0;
394     if (MASTER(cr))
395     {
396         isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
397         if (rerun_fr.natoms != top_global->natoms)
398         {
399             gmx_fatal(FARGS,
400                       "Number of atoms in trajectory (%d) does not match the "
401                       "run input file (%d)\n",
402                       rerun_fr.natoms, top_global->natoms);
403         }
404
405         if (ir->pbcType != PbcType::No)
406         {
407             if (!rerun_fr.bBox)
408             {
409                 gmx_fatal(FARGS,
410                           "Rerun trajectory frame step %" PRId64
411                           " time %f "
412                           "does not contain a box, while pbc is used",
413                           rerun_fr.step, rerun_fr.time);
414             }
415             if (max_cutoff2(ir->pbcType, rerun_fr.box) < gmx::square(fr->rlist))
416             {
417                 gmx_fatal(FARGS,
418                           "Rerun trajectory frame step %" PRId64
419                           " time %f "
420                           "has too small box dimensions",
421                           rerun_fr.step, rerun_fr.time);
422             }
423         }
424     }
425
426     GMX_LOG(mdlog.info)
427             .asParagraph()
428             .appendText(
429                     "Rerun does not report kinetic energy, total energy, temperature, virial and "
430                     "pressure.");
431
432     if (PAR(cr))
433     {
434         rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
435     }
436
437     if (ir->pbcType != PbcType::No)
438     {
439         /* Set the shift vectors.
440          * Necessary here when have a static box different from the tpr box.
441          */
442         calc_shifts(rerun_fr.box, fr->shift_vec);
443     }
444
445     step     = ir->init_step;
446     step_rel = 0;
447
448     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
449             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false, MASTER(cr),
450             ir->nstlist, mdrunOptions.reproducible, nstglobalcomm, mdrunOptions.maximumHoursToRun,
451             ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
452
453     // we don't do counter resetting in rerun - finish will always be valid
454     walltime_accounting_set_valid_finish(walltime_accounting);
455
456     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
457
458     /* and stop now if we should */
459     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
460     while (!isLastStep)
461     {
462         wallcycle_start(wcycle, ewcSTEP);
463
464         if (rerun_fr.bStep)
465         {
466             step     = rerun_fr.step;
467             step_rel = step - ir->init_step;
468         }
469         if (rerun_fr.bTime)
470         {
471             t = rerun_fr.time;
472         }
473         else
474         {
475             t = step;
476         }
477
478         if (ir->efep != efepNO && MASTER(cr))
479         {
480             setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
481         }
482
483         if (MASTER(cr))
484         {
485             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
486             if (constructVsites && DOMAINDECOMP(cr))
487             {
488                 gmx_fatal(FARGS,
489                           "Vsite recalculation with -rerun is not implemented with domain "
490                           "decomposition, "
491                           "use a single rank");
492             }
493             prepareRerunState(rerun_fr, state_global, constructVsites, vsite, ir->delta_t);
494         }
495
496         isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
497
498         if (DOMAINDECOMP(cr))
499         {
500             /* Repartition the domain decomposition */
501             const bool bMasterState = true;
502             dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
503                                 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
504                                 fr, vsite, constr, nrnb, wcycle, mdrunOptions.verbose);
505             shouldCheckNumberOfBondedInteractions = true;
506         }
507
508         if (MASTER(cr))
509         {
510             EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
511         }
512
513         if (ir->efep != efepNO)
514         {
515             update_mdatoms(mdatoms, state->lambda[efptMASS]);
516         }
517
518         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
519                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
520                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
521
522         if (shellfc)
523         {
524             /* Now is the time to relax the shells */
525             relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
526                                 imdSession, pull_work, bNS, force_flags, &top, constr, enerd,
527                                 state->natoms, state->x.arrayRefWithPadding(),
528                                 state->v.arrayRefWithPadding(), state->box, state->lambda,
529                                 &state->hist, &f.view(), force_vir, mdatoms, nrnb, wcycle, shellfc,
530                                 fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
531         }
532         else
533         {
534             /* The coordinates (x) are shifted (to get whole molecules)
535              * in do_force.
536              * This is parallellized as well, and does communication too.
537              * Check comments in sim_util.c
538              */
539             Awh*       awh = nullptr;
540             gmx_edsam* ed  = nullptr;
541             do_force(fplog, cr, ms, ir, awh, enforcedRotation, imdSession, pull_work, step, nrnb,
542                      wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
543                      &f.view(), force_vir, mdatoms, enerd, state->lambda, fr, runScheduleWork,
544                      vsite, mu_tot, t, ed, GMX_FORCE_NS | force_flags, 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.view().force(), isCheckpointingStep,
557                                      doRerun, 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, mdatoms, nrnb, vcm, wcycle,
578                             enerd, force_vir, shake_vir, total_vir, pres, constr, &signaller,
579                             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         /* Output stuff */
594         if (MASTER(cr))
595         {
596             const bool bCalcEnerStep = true;
597             energyOutput.addDataAtEnergyStep(
598                     doFreeEnergyPerturbation, bCalcEnerStep, t, mdatoms->tmass, enerd, ir->fepvals,
599                     ir->expandedvals, state->box,
600                     PTCouplingArrays({ state->boxv, state->nosehoover_xi, state->nosehoover_vxi,
601                                        state->nhpres_xi, state->nhpres_vxi }),
602                     state->fep_state, shake_vir, force_vir, total_vir, pres, ekind, mu_tot, constr);
603
604             const bool do_ene = true;
605             const bool do_log = true;
606             Awh*       awh    = nullptr;
607             const bool do_dr  = ir->nstdisreout != 0;
608             const bool do_or  = ir->nstorireout != 0;
609
610             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
611             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
612                                                do_log ? fplog : nullptr, step, t, fr->fcdata.get(), awh);
613
614             if (do_per_step(step, ir->nstlog))
615             {
616                 if (fflush(fplog) != 0)
617                 {
618                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
619                 }
620             }
621         }
622
623         /* Print the remaining wall clock time for the run */
624         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
625         {
626             if (shellfc)
627             {
628                 fprintf(stderr, "\n");
629             }
630             print_time(stderr, walltime_accounting, step, ir, cr);
631         }
632
633         /* Ion/water position swapping.
634          * Not done in last step since trajectory writing happens before this call
635          * in the MD loop and exchanges would be lost anyway. */
636         if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep && do_per_step(step, ir->swap->nstswap))
637         {
638             const bool doRerun = true;
639             do_swapcoords(cr, step, t, ir, swap, wcycle, rerun_fr.x, rerun_fr.box,
640                           MASTER(cr) && mdrunOptions.verbose, doRerun);
641         }
642
643         if (MASTER(cr))
644         {
645             /* read next frame from input trajectory */
646             isLastStep = !read_next_frame(oenv, status, &rerun_fr);
647         }
648
649         if (PAR(cr))
650         {
651             rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
652         }
653
654         cycles = wallcycle_stop(wcycle, ewcSTEP);
655         if (DOMAINDECOMP(cr) && wcycle)
656         {
657             dd_cycles_add(cr->dd, cycles, ddCyclStep);
658         }
659
660         if (!rerun_fr.bStep)
661         {
662             /* increase the MD step number */
663             step++;
664             step_rel++;
665         }
666     }
667     /* End of main MD loop */
668
669     /* Closing TNG files can include compressing data. Therefore it is good to do that
670      * before stopping the time measurements. */
671     mdoutf_tng_close(outf);
672
673     /* Stop measuring walltime */
674     walltime_accounting_end_time(walltime_accounting);
675
676     if (MASTER(cr))
677     {
678         close_trx(status);
679     }
680
681     if (!thisRankHasDuty(cr, DUTY_PME))
682     {
683         /* Tell the PME only node to finish */
684         gmx_pme_send_finish(cr);
685     }
686
687     done_mdoutf(outf);
688
689     done_shellfc(fplog, shellfc, step_rel);
690
691     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
692 }