407c0190f6daf49aeacbfee850c1ac53e22d805c
[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/localtopologychecker.h"
60 #include "gromacs/domdec/mdsetup.h"
61 #include "gromacs/domdec/partition.h"
62 #include "gromacs/essentialdynamics/edsam.h"
63 #include "gromacs/ewald/pme_load_balancing.h"
64 #include "gromacs/ewald/pme_pp.h"
65 #include "gromacs/fileio/trxio.h"
66 #include "gromacs/gmxlib/network.h"
67 #include "gromacs/gmxlib/nrnb.h"
68 #include "gromacs/gpu_utils/gpu_utils.h"
69 #include "gromacs/listed_forces/listed_forces.h"
70 #include "gromacs/math/functions.h"
71 #include "gromacs/math/utilities.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/math/vectypes.h"
74 #include "gromacs/mdlib/checkpointhandler.h"
75 #include "gromacs/mdlib/compute_io.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/ebin.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/energyoutput.h"
80 #include "gromacs/mdlib/expanded.h"
81 #include "gromacs/mdlib/force.h"
82 #include "gromacs/mdlib/force_flags.h"
83 #include "gromacs/mdlib/forcerec.h"
84 #include "gromacs/mdlib/freeenergyparameters.h"
85 #include "gromacs/mdlib/md_support.h"
86 #include "gromacs/mdlib/mdatoms.h"
87 #include "gromacs/mdlib/mdoutf.h"
88 #include "gromacs/mdlib/membed.h"
89 #include "gromacs/mdlib/resethandler.h"
90 #include "gromacs/mdlib/sighandler.h"
91 #include "gromacs/mdlib/simulationsignal.h"
92 #include "gromacs/mdlib/stat.h"
93 #include "gromacs/mdlib/stophandler.h"
94 #include "gromacs/mdlib/tgroup.h"
95 #include "gromacs/mdlib/trajectory_writing.h"
96 #include "gromacs/mdlib/update.h"
97 #include "gromacs/mdlib/vcm.h"
98 #include "gromacs/mdlib/vsite.h"
99 #include "gromacs/mdrunutility/handlerestart.h"
100 #include "gromacs/mdrunutility/multisim.h"
101 #include "gromacs/mdrunutility/printtime.h"
102 #include "gromacs/mdtypes/awh_history.h"
103 #include "gromacs/mdtypes/awh_params.h"
104 #include "gromacs/mdtypes/commrec.h"
105 #include "gromacs/mdtypes/df_history.h"
106 #include "gromacs/mdtypes/energyhistory.h"
107 #include "gromacs/mdtypes/forcebuffers.h"
108 #include "gromacs/mdtypes/forcerec.h"
109 #include "gromacs/mdtypes/group.h"
110 #include "gromacs/mdtypes/inputrec.h"
111 #include "gromacs/mdtypes/interaction_const.h"
112 #include "gromacs/mdtypes/md_enums.h"
113 #include "gromacs/mdtypes/mdatom.h"
114 #include "gromacs/mdtypes/mdrunoptions.h"
115 #include "gromacs/mdtypes/observableshistory.h"
116 #include "gromacs/mdtypes/simulation_workload.h"
117 #include "gromacs/mdtypes/state.h"
118 #include "gromacs/mimic/utilities.h"
119 #include "gromacs/pbcutil/pbc.h"
120 #include "gromacs/pulling/output.h"
121 #include "gromacs/pulling/pull.h"
122 #include "gromacs/swap/swapcoords.h"
123 #include "gromacs/timing/wallcycle.h"
124 #include "gromacs/timing/walltime_accounting.h"
125 #include "gromacs/topology/atoms.h"
126 #include "gromacs/topology/idef.h"
127 #include "gromacs/topology/mtop_util.h"
128 #include "gromacs/topology/topology.h"
129 #include "gromacs/trajectory/trajectoryframe.h"
130 #include "gromacs/utility/basedefinitions.h"
131 #include "gromacs/utility/cstringutil.h"
132 #include "gromacs/utility/fatalerror.h"
133 #include "gromacs/utility/logger.h"
134 #include "gromacs/utility/real.h"
135
136 #include "legacysimulator.h"
137 #include "replicaexchange.h"
138 #include "shellfc.h"
139
140 using gmx::SimulationSignaller;
141 using gmx::VirtualSitesHandler;
142
143 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
144  *
145  * \param[in]     rerunFrame      The trajectory frame to compute energy/forces for
146  * \param[in,out] globalState     The global state container
147  * \param[in]     constructVsites When true, vsite coordinates are constructed
148  * \param[in]     vsite           Vsite setup, can be nullptr when \p constructVsites = false
149  */
150 static void prepareRerunState(const t_trxframe&          rerunFrame,
151                               t_state*                   globalState,
152                               bool                       constructVsites,
153                               const VirtualSitesHandler* vsite)
154 {
155     auto x      = makeArrayRef(globalState->x);
156     auto rerunX = arrayRefFromArray(reinterpret_cast<gmx::RVec*>(rerunFrame.x), globalState->natoms);
157     std::copy(rerunX.begin(), rerunX.end(), x.begin());
158     copy_mat(rerunFrame.box, globalState->box);
159
160     if (constructVsites)
161     {
162         GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
163
164         vsite->construct(globalState->x, globalState->v, globalState->box, gmx::VSiteOperation::PositionsAndVelocities);
165     }
166 }
167
168 void gmx::LegacySimulator::do_rerun()
169 {
170     // TODO Historically, the EM and MD "integrators" used different
171     // names for the t_inputrec *parameter, but these must have the
172     // same name, now that it's a member of a struct. We use this ir
173     // alias to avoid a large ripple of nearly useless changes.
174     // t_inputrec is being replaced by IMdpOptionsProvider, so this
175     // will go away eventually.
176     const t_inputrec* ir = inputrec;
177     int64_t           step, step_rel;
178     double            t;
179     bool              isLastStep               = false;
180     bool              doFreeEnergyPerturbation = false;
181     unsigned int      force_flags;
182     tensor            force_vir, shake_vir, total_vir, pres;
183     t_trxstatus*      status = nullptr;
184     rvec              mu_tot;
185     t_trxframe        rerun_fr;
186     gmx_localtop_t    top(top_global.ffparams);
187     ForceBuffers      f;
188     gmx_global_stat_t gstat;
189     gmx_shellfc_t*    shellfc;
190
191     double cycles;
192
193     SimulationSignals signals;
194     // Most global communnication stages don't propagate mdrun
195     // signals, and will use this object to achieve that.
196     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
197
198     GMX_LOG(mdlog.info)
199             .asParagraph()
200             .appendText(
201                     "Note that it is planned that the command gmx mdrun -rerun will "
202                     "be available in a different form in a future version of GROMACS, "
203                     "e.g. gmx rerun -f.");
204
205     if (ir->efep != FreeEnergyPerturbationType::No
206         && (mdAtoms->mdatoms()->nMassPerturbed > 0 || (constr && constr->havePerturbedConstraints())))
207     {
208         gmx_fatal(FARGS,
209                   "Perturbed masses or constraints are not supported by rerun. "
210                   "Either make a .tpr without mass and constraint perturbation, "
211                   "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
212     }
213     if (ir->bExpanded)
214     {
215         gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
216     }
217     if (ir->bSimTemp)
218     {
219         gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
220     }
221     if (ir->bDoAwh)
222     {
223         gmx_fatal(FARGS, "AWH not supported by rerun.");
224     }
225     if (replExParams.exchangeInterval > 0)
226     {
227         gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
228     }
229     if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
230     {
231         gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
232     }
233     if (ir->bIMD)
234     {
235         gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
236     }
237     if (isMultiSim(ms))
238     {
239         gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
240     }
241     if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc, [](SimulatedAnnealing i) {
242             return i != SimulatedAnnealing::No;
243         }))
244     {
245         gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
246     }
247
248     /* Rerun can't work if an output file name is the same as the input file name.
249      * If this is the case, the user will get an error telling them what the issue is.
250      */
251     if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0
252         || strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
253     {
254         gmx_fatal(FARGS,
255                   "When using mdrun -rerun, the name of the input trajectory file "
256                   "%s cannot be identical to the name of an output file (whether "
257                   "given explicitly with -o or -x, or by default)",
258                   opt2fn("-rerun", nfile, fnm));
259     }
260
261     /* Settings for rerun */
262     {
263         // TODO: Avoid changing inputrec (#3854)
264         auto* nonConstInputrec               = const_cast<t_inputrec*>(inputrec);
265         nonConstInputrec->nstlist            = 1;
266         nonConstInputrec->nstcalcenergy      = 1;
267         nonConstInputrec->nstxout_compressed = 0;
268     }
269     int        nstglobalcomm = 1;
270     const bool bNS           = true;
271
272     const SimulationGroups* groups = &top_global.groups;
273     if (ir->eI == IntegrationAlgorithm::Mimic)
274     {
275         auto* nonConstGlobalTopology                         = const_cast<gmx_mtop_t*>(&top_global);
276         nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(top_global);
277     }
278     int*                fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
279     gmx::ArrayRef<real> lambda    = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
280     initialize_lambdas(fplog,
281                        ir->efep,
282                        ir->bSimTemp,
283                        *ir->fepvals,
284                        ir->simtempvals->temperatures,
285                        gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
286                        MASTER(cr),
287                        fep_state,
288                        lambda);
289     const bool        simulationsShareState = false;
290     gmx_mdoutf*       outf                  = init_mdoutf(fplog,
291                                    nfile,
292                                    fnm,
293                                    mdrunOptions,
294                                    cr,
295                                    outputProvider,
296                                    mdModulesNotifiers,
297                                    ir,
298                                    top_global,
299                                    oenv,
300                                    wcycle,
301                                    StartingBehavior::NewSimulation,
302                                    simulationsShareState,
303                                    ms);
304     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
305                                    top_global,
306                                    *ir,
307                                    pull_work,
308                                    mdoutf_get_fp_dhdl(outf),
309                                    true,
310                                    StartingBehavior::NewSimulation,
311                                    simulationsShareState,
312                                    mdModulesNotifiers);
313
314     gstat = global_stat_init(ir);
315
316     /* Check for polarizable models and flexible constraints */
317     shellfc = init_shell_flexcon(fplog,
318                                  top_global,
319                                  constr ? constr->numFlexibleConstraints() : 0,
320                                  ir->nstcalcenergy,
321                                  DOMAINDECOMP(cr),
322                                  runScheduleWork->simulationWork.useGpuPme);
323
324     {
325         double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
326         if ((io > 2000) && MASTER(cr))
327         {
328             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
329         }
330     }
331
332     // Local state only becomes valid now.
333     std::unique_ptr<t_state> stateInstance;
334     t_state*                 state;
335
336     if (DOMAINDECOMP(cr))
337     {
338         stateInstance = std::make_unique<t_state>();
339         state         = stateInstance.get();
340         dd_init_local_state(*cr->dd, state_global, state);
341
342         /* Distribute the charge groups over the nodes from the master node */
343         dd_partition_system(fplog,
344                             mdlog,
345                             ir->init_step,
346                             cr,
347                             TRUE,
348                             1,
349                             state_global,
350                             top_global,
351                             *ir,
352                             imdSession,
353                             pull_work,
354                             state,
355                             &f,
356                             mdAtoms,
357                             &top,
358                             fr,
359                             vsite,
360                             constr,
361                             nrnb,
362                             nullptr,
363                             FALSE);
364     }
365     else
366     {
367         state_change_natoms(state_global, state_global->natoms);
368         /* Copy the pointer to the global state */
369         state = state_global;
370
371         mdAlgorithmsSetupAtomData(cr, *ir, top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
372     }
373
374     auto* mdatoms = mdAtoms->mdatoms();
375
376     // NOTE: The global state is no longer used at this point.
377     // But state_global is still used as temporary storage space for writing
378     // the global state to file and potentially for replica exchange.
379     // (Global topology should persist.)
380
381     update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
382
383     if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
384     {
385         doFreeEnergyPerturbation = true;
386     }
387
388     {
389         int cglo_flags = CGLO_GSTAT;
390         if (DOMAINDECOMP(cr) && dd_localTopologyChecker(*cr->dd).shouldCheckNumberOfBondedInteractions())
391         {
392             cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
393         }
394         bool   bSumEkinhOld = false;
395         t_vcm* vcm          = nullptr;
396         compute_globals(gstat,
397                         cr,
398                         ir,
399                         fr,
400                         ekind,
401                         makeConstArrayRef(state->x),
402                         makeConstArrayRef(state->v),
403                         state->box,
404                         mdatoms,
405                         nrnb,
406                         vcm,
407                         nullptr,
408                         enerd,
409                         force_vir,
410                         shake_vir,
411                         total_vir,
412                         pres,
413                         gmx::ArrayRef<real>{},
414                         &nullSignaller,
415                         state->box,
416                         &bSumEkinhOld,
417                         cglo_flags);
418         if (DOMAINDECOMP(cr))
419         {
420             dd_localTopologyChecker(cr->dd)->checkNumberOfBondedInteractions(
421                     &top, makeConstArrayRef(state->x), state->box);
422         }
423     }
424
425     if (MASTER(cr))
426     {
427         fprintf(stderr,
428                 "starting md rerun '%s', reading coordinates from"
429                 " input trajectory '%s'\n\n",
430                 *(top_global.name),
431                 opt2fn("-rerun", nfile, fnm));
432         if (mdrunOptions.verbose)
433         {
434             fprintf(stderr,
435                     "Calculated time to finish depends on nsteps from "
436                     "run input file,\nwhich may not correspond to the time "
437                     "needed to process input trajectory.\n\n");
438         }
439         fprintf(fplog, "\n");
440     }
441
442     walltime_accounting_start_time(walltime_accounting);
443     wallcycle_start(wcycle, WallCycleCounter::Run);
444     print_start(fplog, cr, walltime_accounting, "mdrun");
445
446     /***********************************************************
447      *
448      *             Loop over MD steps
449      *
450      ************************************************************/
451
452     if (constr)
453     {
454         GMX_LOG(mdlog.info)
455                 .asParagraph()
456                 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
457     }
458
459     rerun_fr.natoms = 0;
460     if (MASTER(cr))
461     {
462         isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
463         if (rerun_fr.natoms != top_global.natoms)
464         {
465             gmx_fatal(FARGS,
466                       "Number of atoms in trajectory (%d) does not match the "
467                       "run input file (%d)\n",
468                       rerun_fr.natoms,
469                       top_global.natoms);
470         }
471
472         if (ir->pbcType != PbcType::No)
473         {
474             if (!rerun_fr.bBox)
475             {
476                 gmx_fatal(FARGS,
477                           "Rerun trajectory frame step %" PRId64
478                           " time %f "
479                           "does not contain a box, while pbc is used",
480                           rerun_fr.step,
481                           rerun_fr.time);
482             }
483             if (max_cutoff2(ir->pbcType, rerun_fr.box) < gmx::square(fr->rlist))
484             {
485                 gmx_fatal(FARGS,
486                           "Rerun trajectory frame step %" PRId64
487                           " time %f "
488                           "has too small box dimensions",
489                           rerun_fr.step,
490                           rerun_fr.time);
491             }
492         }
493     }
494
495     GMX_LOG(mdlog.info)
496             .asParagraph()
497             .appendText(
498                     "Rerun does not report kinetic energy, total energy, temperature, virial and "
499                     "pressure.");
500
501     if (PAR(cr))
502     {
503         rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
504     }
505
506     if (ir->pbcType != PbcType::No)
507     {
508         /* Set the shift vectors.
509          * Necessary here when have a static box different from the tpr box.
510          */
511         calc_shifts(rerun_fr.box, fr->shift_vec);
512     }
513
514     step     = ir->init_step;
515     step_rel = 0;
516
517     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
518             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
519             false,
520             MASTER(cr),
521             ir->nstlist,
522             mdrunOptions.reproducible,
523             nstglobalcomm,
524             mdrunOptions.maximumHoursToRun,
525             ir->nstlist == 0,
526             fplog,
527             step,
528             bNS,
529             walltime_accounting);
530
531     // we don't do counter resetting in rerun - finish will always be valid
532     walltime_accounting_set_valid_finish(walltime_accounting);
533
534     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
535
536     /* and stop now if we should */
537     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
538     while (!isLastStep)
539     {
540         wallcycle_start(wcycle, WallCycleCounter::Step);
541
542         if (rerun_fr.bStep)
543         {
544             step     = rerun_fr.step;
545             step_rel = step - ir->init_step;
546         }
547         if (rerun_fr.bTime)
548         {
549             t = rerun_fr.time;
550         }
551         else
552         {
553             t = step;
554         }
555
556         if (ir->efep != FreeEnergyPerturbationType::No && MASTER(cr))
557         {
558             if (rerun_fr.bLambda)
559             {
560                 ir->fepvals->init_lambda = rerun_fr.lambda;
561             }
562             else
563             {
564                 if (rerun_fr.bFepState)
565                 {
566                     state->fep_state = rerun_fr.fep_state;
567                 }
568             }
569
570             state_global->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
571         }
572
573         if (MASTER(cr))
574         {
575             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
576             if (constructVsites && DOMAINDECOMP(cr))
577             {
578                 gmx_fatal(FARGS,
579                           "Vsite recalculation with -rerun is not implemented with domain "
580                           "decomposition, "
581                           "use a single rank");
582             }
583             prepareRerunState(rerun_fr, state_global, constructVsites, vsite);
584         }
585
586         isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
587
588         if (DOMAINDECOMP(cr))
589         {
590             /* Repartition the domain decomposition */
591             const bool bMasterState = true;
592             dd_partition_system(fplog,
593                                 mdlog,
594                                 step,
595                                 cr,
596                                 bMasterState,
597                                 nstglobalcomm,
598                                 state_global,
599                                 top_global,
600                                 *ir,
601                                 imdSession,
602                                 pull_work,
603                                 state,
604                                 &f,
605                                 mdAtoms,
606                                 &top,
607                                 fr,
608                                 vsite,
609                                 constr,
610                                 nrnb,
611                                 wcycle,
612                                 mdrunOptions.verbose);
613         }
614
615         if (MASTER(cr))
616         {
617             EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
618         }
619
620         if (ir->efep != FreeEnergyPerturbationType::No)
621         {
622             update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
623         }
624
625         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
626                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
627                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
628
629         if (shellfc)
630         {
631             /* Now is the time to relax the shells */
632             relax_shell_flexcon(fplog,
633                                 cr,
634                                 ms,
635                                 mdrunOptions.verbose,
636                                 enforcedRotation,
637                                 step,
638                                 ir,
639                                 imdSession,
640                                 pull_work,
641                                 bNS,
642                                 force_flags,
643                                 &top,
644                                 constr,
645                                 enerd,
646                                 state->natoms,
647                                 state->x.arrayRefWithPadding(),
648                                 state->v.arrayRefWithPadding(),
649                                 state->box,
650                                 state->lambda,
651                                 &state->hist,
652                                 &f.view(),
653                                 force_vir,
654                                 *mdatoms,
655                                 nrnb,
656                                 wcycle,
657                                 shellfc,
658                                 fr,
659                                 runScheduleWork,
660                                 t,
661                                 mu_tot,
662                                 vsite,
663                                 ddBalanceRegionHandler);
664         }
665         else
666         {
667             /* The coordinates (x) are shifted (to get whole molecules)
668              * in do_force.
669              * This is parallellized as well, and does communication too.
670              * Check comments in sim_util.c
671              */
672             Awh*       awh = nullptr;
673             gmx_edsam* ed  = nullptr;
674             do_force(fplog,
675                      cr,
676                      ms,
677                      *ir,
678                      awh,
679                      enforcedRotation,
680                      imdSession,
681                      pull_work,
682                      step,
683                      nrnb,
684                      wcycle,
685                      &top,
686                      state->box,
687                      state->x.arrayRefWithPadding(),
688                      &state->hist,
689                      &f.view(),
690                      force_vir,
691                      mdatoms,
692                      enerd,
693                      state->lambda,
694                      fr,
695                      runScheduleWork,
696                      vsite,
697                      mu_tot,
698                      t,
699                      ed,
700                      GMX_FORCE_NS | force_flags,
701                      ddBalanceRegionHandler);
702         }
703
704         /* Now we have the energies and forces corresponding to the
705          * coordinates at time t.
706          */
707         {
708             const bool isCheckpointingStep = false;
709             const bool doRerun             = true;
710             const bool bSumEkinhOld        = false;
711             do_md_trajectory_writing(fplog,
712                                      cr,
713                                      nfile,
714                                      fnm,
715                                      step,
716                                      step_rel,
717                                      t,
718                                      ir,
719                                      state,
720                                      state_global,
721                                      observablesHistory,
722                                      top_global,
723                                      fr,
724                                      outf,
725                                      energyOutput,
726                                      ekind,
727                                      f.view().force(),
728                                      isCheckpointingStep,
729                                      doRerun,
730                                      isLastStep,
731                                      mdrunOptions.writeConfout,
732                                      bSumEkinhOld);
733         }
734
735         stopHandler->setSignal();
736
737         {
738             const bool          doInterSimSignal = false;
739             const bool          doIntraSimSignal = true;
740             bool                bSumEkinhOld     = false;
741             t_vcm*              vcm              = nullptr;
742             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
743
744             int cglo_flags = CGLO_GSTAT | CGLO_ENERGY;
745             if (DOMAINDECOMP(cr) && dd_localTopologyChecker(*cr->dd).shouldCheckNumberOfBondedInteractions())
746             {
747                 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
748             }
749             compute_globals(gstat,
750                             cr,
751                             ir,
752                             fr,
753                             ekind,
754                             makeConstArrayRef(state->x),
755                             makeConstArrayRef(state->v),
756                             state->box,
757                             mdatoms,
758                             nrnb,
759                             vcm,
760                             wcycle,
761                             enerd,
762                             force_vir,
763                             shake_vir,
764                             total_vir,
765                             pres,
766                             constr != nullptr ? constr->rmsdData() : gmx::ArrayRef<real>{},
767                             &signaller,
768                             state->box,
769                             &bSumEkinhOld,
770                             cglo_flags);
771             if (DOMAINDECOMP(cr))
772             {
773                 dd_localTopologyChecker(cr->dd)->checkNumberOfBondedInteractions(
774                         &top, makeConstArrayRef(state->x), state->box);
775             }
776         }
777
778         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
779            the virial that should probably be addressed eventually. state->veta has better properies,
780            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
781            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
782
783         /* Output stuff */
784         if (MASTER(cr))
785         {
786             const bool bCalcEnerStep = true;
787             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
788                                              bCalcEnerStep,
789                                              t,
790                                              mdatoms->tmass,
791                                              enerd,
792                                              ir->fepvals.get(),
793                                              ir->expandedvals.get(),
794                                              state->box,
795                                              PTCouplingArrays({ state->boxv,
796                                                                 state->nosehoover_xi,
797                                                                 state->nosehoover_vxi,
798                                                                 state->nhpres_xi,
799                                                                 state->nhpres_vxi }),
800                                              state->fep_state,
801                                              total_vir,
802                                              pres,
803                                              ekind,
804                                              mu_tot,
805                                              constr);
806
807             const bool do_ene = true;
808             const bool do_log = true;
809             Awh*       awh    = nullptr;
810             const bool do_dr  = ir->nstdisreout != 0;
811             const bool do_or  = ir->nstorireout != 0;
812
813             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
814             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
815                                                do_ene,
816                                                do_dr,
817                                                do_or,
818                                                do_log ? fplog : nullptr,
819                                                step,
820                                                t,
821                                                fr->fcdata.get(),
822                                                awh);
823
824             if (ir->bPull)
825             {
826                 pull_print_output(pull_work, step, t);
827             }
828
829             if (do_per_step(step, ir->nstlog))
830             {
831                 if (fflush(fplog) != 0)
832                 {
833                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
834                 }
835             }
836         }
837
838         /* Print the remaining wall clock time for the run */
839         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
840         {
841             if (shellfc)
842             {
843                 fprintf(stderr, "\n");
844             }
845             print_time(stderr, walltime_accounting, step, ir, cr);
846         }
847
848         /* Ion/water position swapping.
849          * Not done in last step since trajectory writing happens before this call
850          * in the MD loop and exchanges would be lost anyway. */
851         if ((ir->eSwapCoords != SwapType::No) && (step > 0) && !isLastStep
852             && do_per_step(step, ir->swap->nstswap))
853         {
854             const bool doRerun = true;
855             do_swapcoords(cr,
856                           step,
857                           t,
858                           ir,
859                           swap,
860                           wcycle,
861                           rerun_fr.x,
862                           rerun_fr.box,
863                           MASTER(cr) && mdrunOptions.verbose,
864                           doRerun);
865         }
866
867         if (MASTER(cr))
868         {
869             /* read next frame from input trajectory */
870             isLastStep = !read_next_frame(oenv, status, &rerun_fr);
871         }
872
873         if (PAR(cr))
874         {
875             rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
876         }
877
878         cycles = wallcycle_stop(wcycle, WallCycleCounter::Step);
879         if (DOMAINDECOMP(cr) && wcycle)
880         {
881             dd_cycles_add(cr->dd, cycles, ddCyclStep);
882         }
883
884         if (!rerun_fr.bStep)
885         {
886             /* increase the MD step number */
887             step++;
888             step_rel++;
889         }
890     }
891     /* End of main MD loop */
892
893     /* Closing TNG files can include compressing data. Therefore it is good to do that
894      * before stopping the time measurements. */
895     mdoutf_tng_close(outf);
896
897     /* Stop measuring walltime */
898     walltime_accounting_end_time(walltime_accounting);
899
900     if (MASTER(cr))
901     {
902         close_trx(status);
903     }
904
905     if (!thisRankHasDuty(cr, DUTY_PME))
906     {
907         /* Tell the PME only node to finish */
908         gmx_pme_send_finish(cr);
909     }
910
911     done_mdoutf(outf);
912
913     done_shellfc(fplog, shellfc, step_rel);
914
915     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
916 }