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