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