cd648a504e98e7b9723bc6f5a010d7dc6c25a250
[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     double            t;
178     bool              isLastStep               = false;
179     bool              doFreeEnergyPerturbation = false;
180     unsigned int      force_flags;
181     tensor            force_vir, shake_vir, total_vir, pres;
182     t_trxstatus*      status = nullptr;
183     rvec              mu_tot;
184     t_trxframe        rerun_fr;
185     gmx_localtop_t    top(top_global.ffparams);
186     ForceBuffers      f;
187     gmx_global_stat_t gstat;
188     gmx_shellfc_t*    shellfc;
189
190     double cycles;
191
192     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     const SimulationGroups* groups = &top_global.groups;
272     if (ir->eI == IntegrationAlgorithm::Mimic)
273     {
274         auto* nonConstGlobalTopology                         = const_cast<gmx_mtop_t*>(&top_global);
275         nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(top_global);
276     }
277     int*                fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
278     gmx::ArrayRef<real> lambda    = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
279     initialize_lambdas(fplog,
280                        ir->efep,
281                        ir->bSimTemp,
282                        *ir->fepvals,
283                        ir->simtempvals->temperatures,
284                        gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
285                        MASTER(cr),
286                        fep_state,
287                        lambda);
288     const bool        simulationsShareState = false;
289     gmx_mdoutf*       outf                  = init_mdoutf(fplog,
290                                    nfile,
291                                    fnm,
292                                    mdrunOptions,
293                                    cr,
294                                    outputProvider,
295                                    mdModulesNotifiers,
296                                    ir,
297                                    top_global,
298                                    oenv,
299                                    wcycle,
300                                    StartingBehavior::NewSimulation,
301                                    simulationsShareState,
302                                    ms);
303     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
304                                    top_global,
305                                    *ir,
306                                    pull_work,
307                                    mdoutf_get_fp_dhdl(outf),
308                                    true,
309                                    StartingBehavior::NewSimulation,
310                                    simulationsShareState,
311                                    mdModulesNotifiers);
312
313     gstat = global_stat_init(ir);
314
315     /* Check for polarizable models and flexible constraints */
316     shellfc = init_shell_flexcon(fplog,
317                                  top_global,
318                                  constr ? constr->numFlexibleConstraints() : 0,
319                                  ir->nstcalcenergy,
320                                  DOMAINDECOMP(cr),
321                                  runScheduleWork->simulationWork.useGpuPme);
322
323     {
324         double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
325         if ((io > 2000) && MASTER(cr))
326         {
327             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
328         }
329     }
330
331     // Local state only becomes valid now.
332     std::unique_ptr<t_state> stateInstance;
333     t_state*                 state;
334
335     if (DOMAINDECOMP(cr))
336     {
337         stateInstance = std::make_unique<t_state>();
338         state         = stateInstance.get();
339         dd_init_local_state(*cr->dd, state_global, state);
340
341         /* Distribute the charge groups over the nodes from the master node */
342         dd_partition_system(fplog,
343                             mdlog,
344                             ir->init_step,
345                             cr,
346                             TRUE,
347                             1,
348                             state_global,
349                             top_global,
350                             *ir,
351                             imdSession,
352                             pull_work,
353                             state,
354                             &f,
355                             mdAtoms,
356                             &top,
357                             fr,
358                             vsite,
359                             constr,
360                             nrnb,
361                             nullptr,
362                             FALSE);
363     }
364     else
365     {
366         state_change_natoms(state_global, state_global->natoms);
367         /* Copy the pointer to the global state */
368         state = state_global;
369
370         mdAlgorithmsSetupAtomData(cr, *ir, top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
371     }
372
373     auto* mdatoms = mdAtoms->mdatoms();
374
375     // NOTE: The global state is no longer used at this point.
376     // But state_global is still used as temporary storage space for writing
377     // the global state to file and potentially for replica exchange.
378     // (Global topology should persist.)
379
380     update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
381
382     if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
383     {
384         doFreeEnergyPerturbation = true;
385     }
386
387     {
388         int cglo_flags = CGLO_GSTAT;
389         if (DOMAINDECOMP(cr) && dd_localTopologyChecker(*cr->dd).shouldCheckNumberOfBondedInteractions())
390         {
391             cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
392         }
393         bool   bSumEkinhOld = false;
394         t_vcm* vcm          = nullptr;
395         compute_globals(gstat,
396                         cr,
397                         ir,
398                         fr,
399                         ekind,
400                         makeConstArrayRef(state->x),
401                         makeConstArrayRef(state->v),
402                         state->box,
403                         mdatoms,
404                         nrnb,
405                         vcm,
406                         nullptr,
407                         enerd,
408                         force_vir,
409                         shake_vir,
410                         total_vir,
411                         pres,
412                         gmx::ArrayRef<real>{},
413                         &nullSignaller,
414                         state->box,
415                         &bSumEkinhOld,
416                         cglo_flags);
417         if (DOMAINDECOMP(cr))
418         {
419             dd_localTopologyChecker(cr->dd)->checkNumberOfBondedInteractions(
420                     &top, makeConstArrayRef(state->x), state->box);
421         }
422     }
423
424     if (MASTER(cr))
425     {
426         fprintf(stderr,
427                 "starting md rerun '%s', reading coordinates from"
428                 " input trajectory '%s'\n\n",
429                 *(top_global.name),
430                 opt2fn("-rerun", nfile, fnm));
431         if (mdrunOptions.verbose)
432         {
433             fprintf(stderr,
434                     "Calculated time to finish depends on nsteps from "
435                     "run input file,\nwhich may not correspond to the time "
436                     "needed to process input trajectory.\n\n");
437         }
438         fprintf(fplog, "\n");
439     }
440
441     walltime_accounting_start_time(walltime_accounting);
442     wallcycle_start(wcycle, WallCycleCounter::Run);
443     print_start(fplog, cr, walltime_accounting, "mdrun");
444
445     /***********************************************************
446      *
447      *             Loop over MD steps
448      *
449      ************************************************************/
450
451     if (constr)
452     {
453         GMX_LOG(mdlog.info)
454                 .asParagraph()
455                 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
456     }
457
458     rerun_fr.natoms = 0;
459     if (MASTER(cr))
460     {
461         isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
462         if (rerun_fr.natoms != top_global.natoms)
463         {
464             gmx_fatal(FARGS,
465                       "Number of atoms in trajectory (%d) does not match the "
466                       "run input file (%d)\n",
467                       rerun_fr.natoms,
468                       top_global.natoms);
469         }
470
471         if (ir->pbcType != PbcType::No)
472         {
473             if (!rerun_fr.bBox)
474             {
475                 gmx_fatal(FARGS,
476                           "Rerun trajectory frame step %" PRId64
477                           " time %f "
478                           "does not contain a box, while pbc is used",
479                           rerun_fr.step,
480                           rerun_fr.time);
481             }
482             if (max_cutoff2(ir->pbcType, rerun_fr.box) < gmx::square(fr->rlist))
483             {
484                 gmx_fatal(FARGS,
485                           "Rerun trajectory frame step %" PRId64
486                           " time %f "
487                           "has too small box dimensions",
488                           rerun_fr.step,
489                           rerun_fr.time);
490             }
491         }
492     }
493
494     GMX_LOG(mdlog.info)
495             .asParagraph()
496             .appendText(
497                     "Rerun does not report kinetic energy, total energy, temperature, virial and "
498                     "pressure.");
499
500     if (PAR(cr))
501     {
502         rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
503     }
504
505     if (ir->pbcType != PbcType::No)
506     {
507         /* Set the shift vectors.
508          * Necessary here when have a static box different from the tpr box.
509          */
510         calc_shifts(rerun_fr.box, fr->shift_vec);
511     }
512
513     int64_t step     = ir->init_step;
514     int64_t step_rel = 0;
515
516     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
517             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
518             false,
519             MASTER(cr),
520             ir->nstlist,
521             mdrunOptions.reproducible,
522             nstglobalcomm,
523             mdrunOptions.maximumHoursToRun,
524             ir->nstlist == 0,
525             fplog,
526             step,
527             bNS,
528             walltime_accounting);
529
530     // we don't do counter resetting in rerun - finish will always be valid
531     walltime_accounting_set_valid_finish(walltime_accounting);
532
533     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
534
535     /* and stop now if we should */
536     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
537     while (!isLastStep)
538     {
539         wallcycle_start(wcycle, WallCycleCounter::Step);
540
541         if (rerun_fr.bStep)
542         {
543             step     = rerun_fr.step;
544             step_rel = step - ir->init_step;
545         }
546         if (rerun_fr.bTime)
547         {
548             t = rerun_fr.time;
549         }
550         else
551         {
552             t = step;
553         }
554
555         if (ir->efep != FreeEnergyPerturbationType::No && MASTER(cr))
556         {
557             if (rerun_fr.bLambda)
558             {
559                 ir->fepvals->init_lambda = rerun_fr.lambda;
560             }
561             else
562             {
563                 if (rerun_fr.bFepState)
564                 {
565                     state->fep_state = rerun_fr.fep_state;
566                 }
567             }
568
569             state_global->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
570         }
571
572         if (MASTER(cr))
573         {
574             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
575             if (constructVsites && DOMAINDECOMP(cr))
576             {
577                 gmx_fatal(FARGS,
578                           "Vsite recalculation with -rerun is not implemented with domain "
579                           "decomposition, "
580                           "use a single rank");
581             }
582             prepareRerunState(rerun_fr, state_global, constructVsites, vsite);
583         }
584
585         isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
586
587         if (DOMAINDECOMP(cr))
588         {
589             /* Repartition the domain decomposition */
590             const bool bMasterState = true;
591             dd_partition_system(fplog,
592                                 mdlog,
593                                 step,
594                                 cr,
595                                 bMasterState,
596                                 nstglobalcomm,
597                                 state_global,
598                                 top_global,
599                                 *ir,
600                                 imdSession,
601                                 pull_work,
602                                 state,
603                                 &f,
604                                 mdAtoms,
605                                 &top,
606                                 fr,
607                                 vsite,
608                                 constr,
609                                 nrnb,
610                                 wcycle,
611                                 mdrunOptions.verbose);
612         }
613
614         if (MASTER(cr))
615         {
616             EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
617         }
618
619         if (ir->efep != FreeEnergyPerturbationType::No)
620         {
621             update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
622         }
623
624         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
625                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
626                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
627
628         if (shellfc)
629         {
630             /* Now is the time to relax the shells */
631             relax_shell_flexcon(fplog,
632                                 cr,
633                                 ms,
634                                 mdrunOptions.verbose,
635                                 enforcedRotation,
636                                 step,
637                                 ir,
638                                 imdSession,
639                                 pull_work,
640                                 bNS,
641                                 force_flags,
642                                 &top,
643                                 constr,
644                                 enerd,
645                                 state->natoms,
646                                 state->x.arrayRefWithPadding(),
647                                 state->v.arrayRefWithPadding(),
648                                 state->box,
649                                 state->lambda,
650                                 &state->hist,
651                                 &f.view(),
652                                 force_vir,
653                                 *mdatoms,
654                                 nrnb,
655                                 wcycle,
656                                 shellfc,
657                                 fr,
658                                 runScheduleWork,
659                                 t,
660                                 mu_tot,
661                                 vsite,
662                                 ddBalanceRegionHandler);
663         }
664         else
665         {
666             /* The coordinates (x) are shifted (to get whole molecules)
667              * in do_force.
668              * This is parallellized as well, and does communication too.
669              * Check comments in sim_util.c
670              */
671             Awh*       awh = nullptr;
672             gmx_edsam* ed  = nullptr;
673             do_force(fplog,
674                      cr,
675                      ms,
676                      *ir,
677                      awh,
678                      enforcedRotation,
679                      imdSession,
680                      pull_work,
681                      step,
682                      nrnb,
683                      wcycle,
684                      &top,
685                      state->box,
686                      state->x.arrayRefWithPadding(),
687                      &state->hist,
688                      &f.view(),
689                      force_vir,
690                      mdatoms,
691                      enerd,
692                      state->lambda,
693                      fr,
694                      runScheduleWork,
695                      vsite,
696                      mu_tot,
697                      t,
698                      ed,
699                      GMX_FORCE_NS | force_flags,
700                      ddBalanceRegionHandler);
701         }
702
703         /* Now we have the energies and forces corresponding to the
704          * coordinates at time t.
705          */
706         {
707             const bool isCheckpointingStep = false;
708             const bool doRerun             = true;
709             const bool bSumEkinhOld        = false;
710             do_md_trajectory_writing(fplog,
711                                      cr,
712                                      nfile,
713                                      fnm,
714                                      step,
715                                      step_rel,
716                                      t,
717                                      ir,
718                                      state,
719                                      state_global,
720                                      observablesHistory,
721                                      top_global,
722                                      fr,
723                                      outf,
724                                      energyOutput,
725                                      ekind,
726                                      f.view().force(),
727                                      isCheckpointingStep,
728                                      doRerun,
729                                      isLastStep,
730                                      mdrunOptions.writeConfout,
731                                      bSumEkinhOld);
732         }
733
734         stopHandler->setSignal();
735
736         {
737             const bool          doInterSimSignal = false;
738             const bool          doIntraSimSignal = true;
739             bool                bSumEkinhOld     = false;
740             t_vcm*              vcm              = nullptr;
741             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
742
743             int cglo_flags = CGLO_GSTAT | CGLO_ENERGY;
744             if (DOMAINDECOMP(cr) && dd_localTopologyChecker(*cr->dd).shouldCheckNumberOfBondedInteractions())
745             {
746                 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
747             }
748             compute_globals(gstat,
749                             cr,
750                             ir,
751                             fr,
752                             ekind,
753                             makeConstArrayRef(state->x),
754                             makeConstArrayRef(state->v),
755                             state->box,
756                             mdatoms,
757                             nrnb,
758                             vcm,
759                             wcycle,
760                             enerd,
761                             force_vir,
762                             shake_vir,
763                             total_vir,
764                             pres,
765                             constr != nullptr ? constr->rmsdData() : gmx::ArrayRef<real>{},
766                             &signaller,
767                             state->box,
768                             &bSumEkinhOld,
769                             cglo_flags);
770             if (DOMAINDECOMP(cr))
771             {
772                 dd_localTopologyChecker(cr->dd)->checkNumberOfBondedInteractions(
773                         &top, makeConstArrayRef(state->x), state->box);
774             }
775         }
776
777         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
778            the virial that should probably be addressed eventually. state->veta has better properies,
779            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
780            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
781
782         /* Output stuff */
783         if (MASTER(cr))
784         {
785             const bool bCalcEnerStep = true;
786             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
787                                              bCalcEnerStep,
788                                              t,
789                                              mdatoms->tmass,
790                                              enerd,
791                                              ir->fepvals.get(),
792                                              ir->expandedvals.get(),
793                                              state->box,
794                                              PTCouplingArrays({ state->boxv,
795                                                                 state->nosehoover_xi,
796                                                                 state->nosehoover_vxi,
797                                                                 state->nhpres_xi,
798                                                                 state->nhpres_vxi }),
799                                              state->fep_state,
800                                              total_vir,
801                                              pres,
802                                              ekind,
803                                              mu_tot,
804                                              constr);
805
806             const bool do_ene = true;
807             const bool do_log = true;
808             Awh*       awh    = nullptr;
809             const bool do_dr  = ir->nstdisreout != 0;
810             const bool do_or  = ir->nstorireout != 0;
811
812             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
813             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
814                                                do_ene,
815                                                do_dr,
816                                                do_or,
817                                                do_log ? fplog : nullptr,
818                                                step,
819                                                t,
820                                                fr->fcdata.get(),
821                                                awh);
822
823             if (ir->bPull)
824             {
825                 pull_print_output(pull_work, step, t);
826             }
827
828             if (do_per_step(step, ir->nstlog))
829             {
830                 if (fflush(fplog) != 0)
831                 {
832                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
833                 }
834             }
835         }
836
837         /* Print the remaining wall clock time for the run */
838         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
839         {
840             if (shellfc)
841             {
842                 fprintf(stderr, "\n");
843             }
844             print_time(stderr, walltime_accounting, step, ir, cr);
845         }
846
847         /* Ion/water position swapping.
848          * Not done in last step since trajectory writing happens before this call
849          * in the MD loop and exchanges would be lost anyway. */
850         if ((ir->eSwapCoords != SwapType::No) && (step > 0) && !isLastStep
851             && do_per_step(step, ir->swap->nstswap))
852         {
853             const bool doRerun = true;
854             do_swapcoords(cr,
855                           step,
856                           t,
857                           ir,
858                           swap,
859                           wcycle,
860                           rerun_fr.x,
861                           rerun_fr.box,
862                           MASTER(cr) && mdrunOptions.verbose,
863                           doRerun);
864         }
865
866         if (MASTER(cr))
867         {
868             /* read next frame from input trajectory */
869             isLastStep = !read_next_frame(oenv, status, &rerun_fr);
870         }
871
872         if (PAR(cr))
873         {
874             rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
875         }
876
877         cycles = wallcycle_stop(wcycle, WallCycleCounter::Step);
878         if (DOMAINDECOMP(cr) && wcycle)
879         {
880             dd_cycles_add(cr->dd, cycles, ddCyclStep);
881         }
882
883         if (!rerun_fr.bStep)
884         {
885             /* increase the MD step number */
886             step++;
887             step_rel++;
888         }
889     }
890     /* End of main MD loop */
891
892     /* Closing TNG files can include compressing data. Therefore it is good to do that
893      * before stopping the time measurements. */
894     mdoutf_tng_close(outf);
895
896     /* Stop measuring walltime */
897     walltime_accounting_end_time(walltime_accounting);
898
899     if (MASTER(cr))
900     {
901         close_trx(status);
902     }
903
904     if (!thisRankHasDuty(cr, DUTY_PME))
905     {
906         /* Tell the PME only node to finish */
907         gmx_pme_send_finish(cr);
908     }
909
910     done_mdoutf(outf);
911
912     done_shellfc(fplog, shellfc, step_rel);
913
914     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
915 }