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