Move virtual site calculation
[alexxy/gromacs.git] / src / gromacs / mdrun / rerun.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  *
37  * \brief Implements the loop for simulation reruns
38  *
39  * \author Pascal Merz <pascal.merz@colorado.edu>
40  * \ingroup module_mdrun
41  */
42 #include "gmxpre.h"
43
44 #include <cinttypes>
45 #include <cmath>
46 #include <cstdio>
47 #include <cstdlib>
48
49 #include <algorithm>
50 #include <memory>
51
52 #include "gromacs/applied_forces/awh/awh.h"
53 #include "gromacs/commandline/filenm.h"
54 #include "gromacs/domdec/collect.h"
55 #include "gromacs/domdec/dlbtiming.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/domdec/mdsetup.h"
60 #include "gromacs/domdec/partition.h"
61 #include "gromacs/essentialdynamics/edsam.h"
62 #include "gromacs/ewald/pme_load_balancing.h"
63 #include "gromacs/ewald/pme_pp.h"
64 #include "gromacs/fileio/trxio.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gmxlib/nrnb.h"
67 #include "gromacs/gpu_utils/gpu_utils.h"
68 #include "gromacs/listed_forces/listed_forces.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/utilities.h"
71 #include "gromacs/math/vec.h"
72 #include "gromacs/math/vectypes.h"
73 #include "gromacs/mdlib/checkpointhandler.h"
74 #include "gromacs/mdlib/compute_io.h"
75 #include "gromacs/mdlib/constr.h"
76 #include "gromacs/mdlib/ebin.h"
77 #include "gromacs/mdlib/enerdata_utils.h"
78 #include "gromacs/mdlib/energyoutput.h"
79 #include "gromacs/mdlib/expanded.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/freeenergyparameters.h"
84 #include "gromacs/mdlib/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/mdoutf.h"
87 #include "gromacs/mdlib/membed.h"
88 #include "gromacs/mdlib/resethandler.h"
89 #include "gromacs/mdlib/sighandler.h"
90 #include "gromacs/mdlib/simulationsignal.h"
91 #include "gromacs/mdlib/stat.h"
92 #include "gromacs/mdlib/stophandler.h"
93 #include "gromacs/mdlib/tgroup.h"
94 #include "gromacs/mdlib/trajectory_writing.h"
95 #include "gromacs/mdlib/update.h"
96 #include "gromacs/mdlib/vcm.h"
97 #include "gromacs/mdlib/vsite.h"
98 #include "gromacs/mdrunutility/handlerestart.h"
99 #include "gromacs/mdrunutility/multisim.h"
100 #include "gromacs/mdrunutility/printtime.h"
101 #include "gromacs/mdtypes/awh_history.h"
102 #include "gromacs/mdtypes/awh_params.h"
103 #include "gromacs/mdtypes/commrec.h"
104 #include "gromacs/mdtypes/df_history.h"
105 #include "gromacs/mdtypes/energyhistory.h"
106 #include "gromacs/mdtypes/forcebuffers.h"
107 #include "gromacs/mdtypes/forcerec.h"
108 #include "gromacs/mdtypes/group.h"
109 #include "gromacs/mdtypes/inputrec.h"
110 #include "gromacs/mdtypes/interaction_const.h"
111 #include "gromacs/mdtypes/md_enums.h"
112 #include "gromacs/mdtypes/mdatom.h"
113 #include "gromacs/mdtypes/mdrunoptions.h"
114 #include "gromacs/mdtypes/observableshistory.h"
115 #include "gromacs/mdtypes/simulation_workload.h"
116 #include "gromacs/mdtypes/state.h"
117 #include "gromacs/mimic/utilities.h"
118 #include "gromacs/pbcutil/pbc.h"
119 #include "gromacs/pulling/pull.h"
120 #include "gromacs/swap/swapcoords.h"
121 #include "gromacs/timing/wallcycle.h"
122 #include "gromacs/timing/walltime_accounting.h"
123 #include "gromacs/topology/atoms.h"
124 #include "gromacs/topology/idef.h"
125 #include "gromacs/topology/mtop_util.h"
126 #include "gromacs/topology/topology.h"
127 #include "gromacs/trajectory/trajectoryframe.h"
128 #include "gromacs/utility/basedefinitions.h"
129 #include "gromacs/utility/cstringutil.h"
130 #include "gromacs/utility/fatalerror.h"
131 #include "gromacs/utility/logger.h"
132 #include "gromacs/utility/real.h"
133
134 #include "legacysimulator.h"
135 #include "replicaexchange.h"
136 #include "shellfc.h"
137
138 using gmx::SimulationSignaller;
139 using gmx::VirtualSitesHandler;
140
141 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
142  *
143  * \param[in]     rerunFrame      The trajectory frame to compute energy/forces for
144  * \param[in,out] globalState     The global state container
145  * \param[in]     constructVsites When true, vsite coordinates are constructed
146  * \param[in]     vsite           Vsite setup, can be nullptr when \p constructVsites = false
147  * \param[in]     timeStep        Time step, used for constructing vsites
148  */
149 static void prepareRerunState(const t_trxframe&          rerunFrame,
150                               t_state*                   globalState,
151                               bool                       constructVsites,
152                               const VirtualSitesHandler* vsite,
153                               double                     timeStep)
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, timeStep, globalState->v, globalState->box);
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     /* Domain decomposition could incorrectly miss a bonded
194        interaction, but checking for that requires a global
195        communication stage, which does not otherwise happen in DD
196        code. So we do that alongside the first global energy reduction
197        after a new DD is made. These variables handle whether the
198        check happens, and the result it returns. */
199     bool shouldCheckNumberOfBondedInteractions = false;
200     int  totalNumberOfBondedInteractions       = -1;
201
202     SimulationSignals signals;
203     // Most global communnication stages don't propagate mdrun
204     // signals, and will use this object to achieve that.
205     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
206
207     GMX_LOG(mdlog.info)
208             .asParagraph()
209             .appendText(
210                     "Note that it is planned that the command gmx mdrun -rerun will "
211                     "be available in a different form in a future version of GROMACS, "
212                     "e.g. gmx rerun -f.");
213
214     if (ir->efep != efepNO
215         && (mdAtoms->mdatoms()->nMassPerturbed > 0 || (constr && constr->havePerturbedConstraints())))
216     {
217         gmx_fatal(FARGS,
218                   "Perturbed masses or constraints are not supported by rerun. "
219                   "Either make a .tpr without mass and constraint perturbation, "
220                   "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
221     }
222     if (ir->bExpanded)
223     {
224         gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
225     }
226     if (ir->bSimTemp)
227     {
228         gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
229     }
230     if (ir->bDoAwh)
231     {
232         gmx_fatal(FARGS, "AWH not supported by rerun.");
233     }
234     if (replExParams.exchangeInterval > 0)
235     {
236         gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
237     }
238     if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
239     {
240         gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
241     }
242     if (ir->bIMD)
243     {
244         gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
245     }
246     if (isMultiSim(ms))
247     {
248         gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
249     }
250     if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc, [](int i) {
251             return i != eannNO;
252         }))
253     {
254         gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
255     }
256
257     /* Rerun can't work if an output file name is the same as the input file name.
258      * If this is the case, the user will get an error telling them what the issue is.
259      */
260     if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0
261         || strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
262     {
263         gmx_fatal(FARGS,
264                   "When using mdrun -rerun, the name of the input trajectory file "
265                   "%s cannot be identical to the name of an output file (whether "
266                   "given explicitly with -o or -x, or by default)",
267                   opt2fn("-rerun", nfile, fnm));
268     }
269
270     /* Settings for rerun */
271     {
272         // TODO: Avoid changing inputrec (#3854)
273         auto* nonConstInputrec               = const_cast<t_inputrec*>(inputrec);
274         nonConstInputrec->nstlist            = 1;
275         nonConstInputrec->nstcalcenergy      = 1;
276         nonConstInputrec->nstxout_compressed = 0;
277     }
278     int        nstglobalcomm = 1;
279     const bool bNS           = true;
280
281     const SimulationGroups* groups = &top_global->groups;
282     if (ir->eI == eiMimic)
283     {
284         auto nonConstGlobalTopology                          = const_cast<gmx_mtop_t*>(top_global);
285         nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(*top_global);
286     }
287     int*                fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
288     gmx::ArrayRef<real> lambda    = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
289     initialize_lambdas(fplog, *ir, MASTER(cr), fep_state, lambda);
290     const bool        simulationsShareState = false;
291     gmx_mdoutf*       outf                  = init_mdoutf(fplog,
292                                    nfile,
293                                    fnm,
294                                    mdrunOptions,
295                                    cr,
296                                    outputProvider,
297                                    mdModulesNotifier,
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                                    mdModulesNotifier);
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     // Local state only becomes valid now.
334     std::unique_ptr<t_state> stateInstance;
335     t_state*                 state;
336
337     if (DOMAINDECOMP(cr))
338     {
339         stateInstance = std::make_unique<t_state>();
340         state         = stateInstance.get();
341         dd_init_local_state(cr->dd, state_global, state);
342
343         /* Distribute the charge groups over the nodes from the master node */
344         dd_partition_system(fplog,
345                             mdlog,
346                             ir->init_step,
347                             cr,
348                             TRUE,
349                             1,
350                             state_global,
351                             *top_global,
352                             ir,
353                             imdSession,
354                             pull_work,
355                             state,
356                             &f,
357                             mdAtoms,
358                             &top,
359                             fr,
360                             vsite,
361                             constr,
362                             nrnb,
363                             nullptr,
364                             FALSE);
365         shouldCheckNumberOfBondedInteractions = true;
366     }
367     else
368     {
369         state_change_natoms(state_global, state_global->natoms);
370         /* Copy the pointer to the global state */
371         state = state_global;
372
373         mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
374     }
375
376     auto mdatoms = mdAtoms->mdatoms();
377
378     // NOTE: The global state is no longer used at this point.
379     // But state_global is still used as temporary storage space for writing
380     // the global state to file and potentially for replica exchange.
381     // (Global topology should persist.)
382
383     update_mdatoms(mdatoms, state->lambda[efptMASS]);
384
385     if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
386     {
387         doFreeEnergyPerturbation = true;
388     }
389
390     {
391         int cglo_flags =
392                 (CGLO_GSTAT
393                  | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
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                         &totalNumberOfBondedInteractions,
417                         &bSumEkinhOld,
418                         cglo_flags);
419     }
420     checkNumberOfBondedInteractions(mdlog,
421                                     cr,
422                                     totalNumberOfBondedInteractions,
423                                     top_global,
424                                     &top,
425                                     makeConstArrayRef(state->x),
426                                     state->box,
427                                     &shouldCheckNumberOfBondedInteractions);
428
429     if (MASTER(cr))
430     {
431         fprintf(stderr,
432                 "starting md rerun '%s', reading coordinates from"
433                 " input trajectory '%s'\n\n",
434                 *(top_global->name),
435                 opt2fn("-rerun", nfile, fnm));
436         if (mdrunOptions.verbose)
437         {
438             fprintf(stderr,
439                     "Calculated time to finish depends on nsteps from "
440                     "run input file,\nwhich may not correspond to the time "
441                     "needed to process input trajectory.\n\n");
442         }
443         fprintf(fplog, "\n");
444     }
445
446     walltime_accounting_start_time(walltime_accounting);
447     wallcycle_start(wcycle, ewcRUN);
448     print_start(fplog, cr, walltime_accounting, "mdrun");
449
450     /***********************************************************
451      *
452      *             Loop over MD steps
453      *
454      ************************************************************/
455
456     if (constr)
457     {
458         GMX_LOG(mdlog.info)
459                 .asParagraph()
460                 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
461     }
462
463     rerun_fr.natoms = 0;
464     if (MASTER(cr))
465     {
466         isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
467         if (rerun_fr.natoms != top_global->natoms)
468         {
469             gmx_fatal(FARGS,
470                       "Number of atoms in trajectory (%d) does not match the "
471                       "run input file (%d)\n",
472                       rerun_fr.natoms,
473                       top_global->natoms);
474         }
475
476         if (ir->pbcType != PbcType::No)
477         {
478             if (!rerun_fr.bBox)
479             {
480                 gmx_fatal(FARGS,
481                           "Rerun trajectory frame step %" PRId64
482                           " time %f "
483                           "does not contain a box, while pbc is used",
484                           rerun_fr.step,
485                           rerun_fr.time);
486             }
487             if (max_cutoff2(ir->pbcType, rerun_fr.box) < gmx::square(fr->rlist))
488             {
489                 gmx_fatal(FARGS,
490                           "Rerun trajectory frame step %" PRId64
491                           " time %f "
492                           "has too small box dimensions",
493                           rerun_fr.step,
494                           rerun_fr.time);
495             }
496         }
497     }
498
499     GMX_LOG(mdlog.info)
500             .asParagraph()
501             .appendText(
502                     "Rerun does not report kinetic energy, total energy, temperature, virial and "
503                     "pressure.");
504
505     if (PAR(cr))
506     {
507         rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
508     }
509
510     if (ir->pbcType != PbcType::No)
511     {
512         /* Set the shift vectors.
513          * Necessary here when have a static box different from the tpr box.
514          */
515         calc_shifts(rerun_fr.box, fr->shift_vec);
516     }
517
518     step     = ir->init_step;
519     step_rel = 0;
520
521     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
522             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
523             false,
524             MASTER(cr),
525             ir->nstlist,
526             mdrunOptions.reproducible,
527             nstglobalcomm,
528             mdrunOptions.maximumHoursToRun,
529             ir->nstlist == 0,
530             fplog,
531             step,
532             bNS,
533             walltime_accounting);
534
535     // we don't do counter resetting in rerun - finish will always be valid
536     walltime_accounting_set_valid_finish(walltime_accounting);
537
538     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
539
540     /* and stop now if we should */
541     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
542     while (!isLastStep)
543     {
544         wallcycle_start(wcycle, ewcSTEP);
545
546         if (rerun_fr.bStep)
547         {
548             step     = rerun_fr.step;
549             step_rel = step - ir->init_step;
550         }
551         if (rerun_fr.bTime)
552         {
553             t = rerun_fr.time;
554         }
555         else
556         {
557             t = step;
558         }
559
560         if (ir->efep != efepNO && MASTER(cr))
561         {
562             if (rerun_fr.bLambda)
563             {
564                 ir->fepvals->init_lambda = rerun_fr.lambda;
565             }
566             else
567             {
568                 if (rerun_fr.bFepState)
569                 {
570                     state->fep_state = rerun_fr.fep_state;
571                 }
572             }
573
574             state_global->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
575         }
576
577         if (MASTER(cr))
578         {
579             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
580             if (constructVsites && DOMAINDECOMP(cr))
581             {
582                 gmx_fatal(FARGS,
583                           "Vsite recalculation with -rerun is not implemented with domain "
584                           "decomposition, "
585                           "use a single rank");
586             }
587             prepareRerunState(rerun_fr, state_global, constructVsites, vsite, ir->delta_t);
588         }
589
590         isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
591
592         if (DOMAINDECOMP(cr))
593         {
594             /* Repartition the domain decomposition */
595             const bool bMasterState = true;
596             dd_partition_system(fplog,
597                                 mdlog,
598                                 step,
599                                 cr,
600                                 bMasterState,
601                                 nstglobalcomm,
602                                 state_global,
603                                 *top_global,
604                                 ir,
605                                 imdSession,
606                                 pull_work,
607                                 state,
608                                 &f,
609                                 mdAtoms,
610                                 &top,
611                                 fr,
612                                 vsite,
613                                 constr,
614                                 nrnb,
615                                 wcycle,
616                                 mdrunOptions.verbose);
617             shouldCheckNumberOfBondedInteractions = true;
618         }
619
620         if (MASTER(cr))
621         {
622             EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
623         }
624
625         if (ir->efep != efepNO)
626         {
627             update_mdatoms(mdatoms, state->lambda[efptMASS]);
628         }
629
630         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
631                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
632                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
633
634         if (shellfc)
635         {
636             /* Now is the time to relax the shells */
637             relax_shell_flexcon(fplog,
638                                 cr,
639                                 ms,
640                                 mdrunOptions.verbose,
641                                 enforcedRotation,
642                                 step,
643                                 ir,
644                                 imdSession,
645                                 pull_work,
646                                 bNS,
647                                 force_flags,
648                                 &top,
649                                 constr,
650                                 enerd,
651                                 state->natoms,
652                                 state->x.arrayRefWithPadding(),
653                                 state->v.arrayRefWithPadding(),
654                                 state->box,
655                                 state->lambda,
656                                 &state->hist,
657                                 &f.view(),
658                                 force_vir,
659                                 mdatoms,
660                                 nrnb,
661                                 wcycle,
662                                 shellfc,
663                                 fr,
664                                 runScheduleWork,
665                                 t,
666                                 mu_tot,
667                                 vsite,
668                                 ddBalanceRegionHandler);
669         }
670         else
671         {
672             /* The coordinates (x) are shifted (to get whole molecules)
673              * in do_force.
674              * This is parallellized as well, and does communication too.
675              * Check comments in sim_util.c
676              */
677             Awh*       awh = nullptr;
678             gmx_edsam* ed  = nullptr;
679             do_force(fplog,
680                      cr,
681                      ms,
682                      ir,
683                      awh,
684                      enforcedRotation,
685                      imdSession,
686                      pull_work,
687                      step,
688                      nrnb,
689                      wcycle,
690                      &top,
691                      state->box,
692                      state->x.arrayRefWithPadding(),
693                      &state->hist,
694                      &f.view(),
695                      force_vir,
696                      mdatoms,
697                      enerd,
698                      state->lambda,
699                      fr,
700                      runScheduleWork,
701                      vsite,
702                      mu_tot,
703                      t,
704                      ed,
705                      GMX_FORCE_NS | force_flags,
706                      ddBalanceRegionHandler);
707         }
708
709         /* Now we have the energies and forces corresponding to the
710          * coordinates at time t.
711          */
712         {
713             const bool isCheckpointingStep = false;
714             const bool doRerun             = true;
715             const bool bSumEkinhOld        = false;
716             do_md_trajectory_writing(fplog,
717                                      cr,
718                                      nfile,
719                                      fnm,
720                                      step,
721                                      step_rel,
722                                      t,
723                                      ir,
724                                      state,
725                                      state_global,
726                                      observablesHistory,
727                                      top_global,
728                                      fr,
729                                      outf,
730                                      energyOutput,
731                                      ekind,
732                                      f.view().force(),
733                                      isCheckpointingStep,
734                                      doRerun,
735                                      isLastStep,
736                                      mdrunOptions.writeConfout,
737                                      bSumEkinhOld);
738         }
739
740         stopHandler->setSignal();
741
742         {
743             const bool          doInterSimSignal = false;
744             const bool          doIntraSimSignal = true;
745             bool                bSumEkinhOld     = false;
746             t_vcm*              vcm              = nullptr;
747             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
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                             &totalNumberOfBondedInteractions,
770                             &bSumEkinhOld,
771                             CGLO_GSTAT | CGLO_ENERGY
772                                     | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
773                                                                              : 0));
774             checkNumberOfBondedInteractions(mdlog,
775                                             cr,
776                                             totalNumberOfBondedInteractions,
777                                             top_global,
778                                             &top,
779                                             makeConstArrayRef(state->x),
780                                             state->box,
781                                             &shouldCheckNumberOfBondedInteractions);
782         }
783
784         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
785            the virial that should probably be addressed eventually. state->veta has better properies,
786            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
787            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
788
789         /* Output stuff */
790         if (MASTER(cr))
791         {
792             const bool bCalcEnerStep = true;
793             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
794                                              bCalcEnerStep,
795                                              t,
796                                              mdatoms->tmass,
797                                              enerd,
798                                              ir->fepvals,
799                                              ir->expandedvals,
800                                              state->box,
801                                              PTCouplingArrays({ state->boxv,
802                                                                 state->nosehoover_xi,
803                                                                 state->nosehoover_vxi,
804                                                                 state->nhpres_xi,
805                                                                 state->nhpres_vxi }),
806                                              state->fep_state,
807                                              shake_vir,
808                                              force_vir,
809                                              total_vir,
810                                              pres,
811                                              ekind,
812                                              mu_tot,
813                                              constr);
814
815             const bool do_ene = true;
816             const bool do_log = true;
817             Awh*       awh    = nullptr;
818             const bool do_dr  = ir->nstdisreout != 0;
819             const bool do_or  = ir->nstorireout != 0;
820
821             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
822             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
823                                                do_ene,
824                                                do_dr,
825                                                do_or,
826                                                do_log ? fplog : nullptr,
827                                                step,
828                                                t,
829                                                fr->fcdata.get(),
830                                                awh);
831
832             if (do_per_step(step, ir->nstlog))
833             {
834                 if (fflush(fplog) != 0)
835                 {
836                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
837                 }
838             }
839         }
840
841         /* Print the remaining wall clock time for the run */
842         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
843         {
844             if (shellfc)
845             {
846                 fprintf(stderr, "\n");
847             }
848             print_time(stderr, walltime_accounting, step, ir, cr);
849         }
850
851         /* Ion/water position swapping.
852          * Not done in last step since trajectory writing happens before this call
853          * in the MD loop and exchanges would be lost anyway. */
854         if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep && do_per_step(step, ir->swap->nstswap))
855         {
856             const bool doRerun = true;
857             do_swapcoords(cr,
858                           step,
859                           t,
860                           ir,
861                           swap,
862                           wcycle,
863                           rerun_fr.x,
864                           rerun_fr.box,
865                           MASTER(cr) && mdrunOptions.verbose,
866                           doRerun);
867         }
868
869         if (MASTER(cr))
870         {
871             /* read next frame from input trajectory */
872             isLastStep = !read_next_frame(oenv, status, &rerun_fr);
873         }
874
875         if (PAR(cr))
876         {
877             rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
878         }
879
880         cycles = wallcycle_stop(wcycle, ewcSTEP);
881         if (DOMAINDECOMP(cr) && wcycle)
882         {
883             dd_cycles_add(cr->dd, cycles, ddCyclStep);
884         }
885
886         if (!rerun_fr.bStep)
887         {
888             /* increase the MD step number */
889             step++;
890             step_rel++;
891         }
892     }
893     /* End of main MD loop */
894
895     /* Closing TNG files can include compressing data. Therefore it is good to do that
896      * before stopping the time measurements. */
897     mdoutf_tng_close(outf);
898
899     /* Stop measuring walltime */
900     walltime_accounting_end_time(walltime_accounting);
901
902     if (MASTER(cr))
903     {
904         close_trx(status);
905     }
906
907     if (!thisRankHasDuty(cr, DUTY_PME))
908     {
909         /* Tell the PME only node to finish */
910         gmx_pme_send_finish(cr);
911     }
912
913     done_mdoutf(outf);
914
915     done_shellfc(fplog, shellfc, step_rel);
916
917     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
918 }