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