ce23edc526ff30eb75d9c175c38209dd5ca70f1a
[alexxy/gromacs.git] / src / gromacs / mdrun / mimic.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
36 /*! \internal \file
37  *
38  * \brief Declares the loop for MiMiC QM/MM
39  *
40  * \author Viacheslav Bolnykh <v.bolnykh@hpc-leap.eu>
41  * \ingroup module_mdrun
42  */
43 #include "gmxpre.h"
44
45 #include <cinttypes>
46 #include <cmath>
47 #include <cstdio>
48 #include <cstdlib>
49
50 #include <algorithm>
51 #include <memory>
52
53 #include "gromacs/applied_forces/awh/awh.h"
54 #include "gromacs/commandline/filenm.h"
55 #include "gromacs/domdec/collect.h"
56 #include "gromacs/domdec/dlbtiming.h"
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_network.h"
59 #include "gromacs/domdec/domdec_struct.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/enerdata.h"
107 #include "gromacs/mdtypes/energyhistory.h"
108 #include "gromacs/mdtypes/forcebuffers.h"
109 #include "gromacs/mdtypes/forcerec.h"
110 #include "gromacs/mdtypes/group.h"
111 #include "gromacs/mdtypes/inputrec.h"
112 #include "gromacs/mdtypes/interaction_const.h"
113 #include "gromacs/mdtypes/md_enums.h"
114 #include "gromacs/mdtypes/mdatom.h"
115 #include "gromacs/mdtypes/mdrunoptions.h"
116 #include "gromacs/mdtypes/observableshistory.h"
117 #include "gromacs/mdtypes/simulation_workload.h"
118 #include "gromacs/mdtypes/state.h"
119 #include "gromacs/mimic/communicator.h"
120 #include "gromacs/mimic/utilities.h"
121 #include "gromacs/pbcutil/pbc.h"
122 #include "gromacs/pulling/pull.h"
123 #include "gromacs/timing/wallcycle.h"
124 #include "gromacs/timing/walltime_accounting.h"
125 #include "gromacs/topology/atoms.h"
126 #include "gromacs/topology/idef.h"
127 #include "gromacs/topology/mtop_util.h"
128 #include "gromacs/topology/topology.h"
129 #include "gromacs/trajectory/trajectoryframe.h"
130 #include "gromacs/utility/basedefinitions.h"
131 #include "gromacs/utility/cstringutil.h"
132 #include "gromacs/utility/fatalerror.h"
133 #include "gromacs/utility/logger.h"
134 #include "gromacs/utility/real.h"
135
136 #include "legacysimulator.h"
137 #include "replicaexchange.h"
138 #include "shellfc.h"
139
140 using gmx::SimulationSignaller;
141
142 void gmx::LegacySimulator::do_mimic()
143 {
144     const t_inputrec* ir = inputrec;
145     int64_t           step, step_rel;
146     double            t;
147     bool              isLastStep               = false;
148     bool              doFreeEnergyPerturbation = false;
149     unsigned int      force_flags;
150     tensor            force_vir, shake_vir, total_vir, pres;
151     rvec              mu_tot;
152     ForceBuffers      f;
153     gmx_global_stat_t gstat;
154     gmx_shellfc_t*    shellfc;
155
156     double cycles;
157
158     SimulationSignals signals;
159     // Most global communnication stages don't propagate mdrun
160     // signals, and will use this object to achieve that.
161     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
162
163     if (ir->bExpanded)
164     {
165         gmx_fatal(FARGS, "Expanded ensemble not supported by MiMiC.");
166     }
167     if (ir->bSimTemp)
168     {
169         gmx_fatal(FARGS, "Simulated tempering not supported by MiMiC.");
170     }
171     if (ir->bDoAwh)
172     {
173         gmx_fatal(FARGS, "AWH not supported by MiMiC.");
174     }
175     if (replExParams.exchangeInterval > 0)
176     {
177         gmx_fatal(FARGS, "Replica exchange not supported by MiMiC.");
178     }
179     if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
180     {
181         gmx_fatal(FARGS, "Essential dynamics not supported by MiMiC.");
182     }
183     if (ir->bIMD)
184     {
185         gmx_fatal(FARGS, "Interactive MD not supported by MiMiC.");
186     }
187     if (isMultiSim(ms))
188     {
189         gmx_fatal(FARGS, "Multiple simulations not supported by MiMiC.");
190     }
191     if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc, [](SimulatedAnnealing i) {
192             return i != SimulatedAnnealing::No;
193         }))
194     {
195         gmx_fatal(FARGS, "Simulated annealing not supported by MiMiC.");
196     }
197
198     /* Settings for rerun */
199     {
200         // TODO: Avoid changing inputrec (#3854)
201         auto* nonConstInputrec               = const_cast<t_inputrec*>(inputrec);
202         nonConstInputrec->nstlist            = 1;
203         nonConstInputrec->nstcalcenergy      = 1;
204         nonConstInputrec->nstxout_compressed = 0;
205     }
206     int        nstglobalcomm = 1;
207     const bool bNS           = true;
208
209     if (MASTER(cr))
210     {
211         MimicCommunicator::init();
212         auto nonConstGlobalTopology = const_cast<gmx_mtop_t*>(&top_global);
213         MimicCommunicator::sendInitData(nonConstGlobalTopology, state_global->x);
214         // TODO: Avoid changing inputrec (#3854)
215         auto* nonConstInputrec   = const_cast<t_inputrec*>(inputrec);
216         nonConstInputrec->nsteps = MimicCommunicator::getStepNumber();
217     }
218     if (DOMAINDECOMP(cr))
219     {
220         // TODO: Avoid changing inputrec (#3854)
221         auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
222         gmx_bcast(sizeof(ir->nsteps), &nonConstInputrec->nsteps, cr->mpi_comm_mygroup);
223     }
224
225     const SimulationGroups* groups = &top_global.groups;
226     {
227         auto nonConstGlobalTopology                          = const_cast<gmx_mtop_t*>(&top_global);
228         nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(top_global);
229     }
230
231     initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda);
232
233     const bool        simulationsShareState = false;
234     gmx_mdoutf*       outf                  = init_mdoutf(fplog,
235                                    nfile,
236                                    fnm,
237                                    mdrunOptions,
238                                    cr,
239                                    outputProvider,
240                                    mdModulesNotifier,
241                                    ir,
242                                    top_global,
243                                    oenv,
244                                    wcycle,
245                                    StartingBehavior::NewSimulation,
246                                    simulationsShareState,
247                                    ms);
248     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
249                                    top_global,
250                                    *ir,
251                                    pull_work,
252                                    mdoutf_get_fp_dhdl(outf),
253                                    true,
254                                    StartingBehavior::NewSimulation,
255                                    simulationsShareState,
256                                    mdModulesNotifier);
257
258     gstat = global_stat_init(ir);
259
260     /* Check for polarizable models and flexible constraints */
261     shellfc = init_shell_flexcon(fplog,
262                                  top_global,
263                                  constr ? constr->numFlexibleConstraints() : 0,
264                                  ir->nstcalcenergy,
265                                  DOMAINDECOMP(cr),
266                                  runScheduleWork->simulationWork.useGpuPme);
267
268     {
269         double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
270         if ((io > 2000) && MASTER(cr))
271         {
272             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
273         }
274     }
275
276     // Local state only becomes valid now.
277     std::unique_ptr<t_state> stateInstance;
278     t_state*                 state;
279
280     gmx_localtop_t top(top_global.ffparams);
281
282     if (DOMAINDECOMP(cr))
283     {
284         stateInstance = std::make_unique<t_state>();
285         state         = stateInstance.get();
286         dd_init_local_state(*cr->dd, state_global, state);
287
288         /* Distribute the charge groups over the nodes from the master node */
289         dd_partition_system(fplog,
290                             mdlog,
291                             ir->init_step,
292                             cr,
293                             TRUE,
294                             1,
295                             state_global,
296                             top_global,
297                             *ir,
298                             imdSession,
299                             pull_work,
300                             state,
301                             &f,
302                             mdAtoms,
303                             &top,
304                             fr,
305                             vsite,
306                             constr,
307                             nrnb,
308                             nullptr,
309                             FALSE);
310     }
311     else
312     {
313         state_change_natoms(state_global, state_global->natoms);
314         /* Copy the pointer to the global state */
315         state = state_global;
316
317         mdAlgorithmsSetupAtomData(cr, *ir, top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
318     }
319
320     auto mdatoms = mdAtoms->mdatoms();
321
322     // NOTE: The global state is no longer used at this point.
323     // But state_global is still used as temporary storage space for writing
324     // the global state to file and potentially for replica exchange.
325     // (Global topology should persist.)
326
327     update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
328
329     if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
330     {
331         doFreeEnergyPerturbation = true;
332     }
333
334     {
335         int cglo_flags = CGLO_GSTAT;
336         if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
337         {
338             cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
339         }
340         bool   bSumEkinhOld = false;
341         t_vcm* vcm          = nullptr;
342         compute_globals(gstat,
343                         cr,
344                         ir,
345                         fr,
346                         ekind,
347                         makeConstArrayRef(state->x),
348                         makeConstArrayRef(state->v),
349                         state->box,
350                         mdatoms,
351                         nrnb,
352                         vcm,
353                         nullptr,
354                         enerd,
355                         force_vir,
356                         shake_vir,
357                         total_vir,
358                         pres,
359                         gmx::ArrayRef<real>{},
360                         &nullSignaller,
361                         state->box,
362                         &bSumEkinhOld,
363                         cglo_flags);
364         if (DOMAINDECOMP(cr))
365         {
366             checkNumberOfBondedInteractions(
367                     mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
368         }
369     }
370
371     if (MASTER(cr))
372     {
373         fprintf(stderr, "starting MiMiC MD run '%s'\n\n", *(top_global.name));
374         if (mdrunOptions.verbose)
375         {
376             fprintf(stderr,
377                     "Calculated time to finish depends on nsteps from "
378                     "run input file,\nwhich may not correspond to the time "
379                     "needed to process input trajectory.\n\n");
380         }
381         fprintf(fplog, "\n");
382     }
383
384     walltime_accounting_start_time(walltime_accounting);
385     wallcycle_start(wcycle, ewcRUN);
386     print_start(fplog, cr, walltime_accounting, "mdrun");
387
388     /***********************************************************
389      *
390      *             Loop over MD steps
391      *
392      ************************************************************/
393
394     if (constr)
395     {
396         GMX_LOG(mdlog.info)
397                 .asParagraph()
398                 .appendText(
399                         "Simulations has constraints. Constraints will "
400                         "be handled by CPMD.");
401     }
402
403     GMX_LOG(mdlog.info)
404             .asParagraph()
405             .appendText(
406                     "MiMiC does not report kinetic energy, total energy, temperature, virial and "
407                     "pressure.");
408
409     step     = ir->init_step;
410     step_rel = 0;
411
412     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
413             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
414             false,
415             MASTER(cr),
416             ir->nstlist,
417             mdrunOptions.reproducible,
418             nstglobalcomm,
419             mdrunOptions.maximumHoursToRun,
420             ir->nstlist == 0,
421             fplog,
422             step,
423             bNS,
424             walltime_accounting);
425
426     // we don't do counter resetting in rerun - finish will always be valid
427     walltime_accounting_set_valid_finish(walltime_accounting);
428
429     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
430
431     /* and stop now if we should */
432     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
433     while (!isLastStep)
434     {
435         isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel == ir->nsteps));
436         wallcycle_start(wcycle, ewcSTEP);
437
438         t = step;
439
440         if (MASTER(cr))
441         {
442             MimicCommunicator::getCoords(&state_global->x, state_global->natoms);
443         }
444
445         if (ir->efep != FreeEnergyPerturbationType::No)
446         {
447             state->lambda = currentLambdas(step, *(ir->fepvals), state_global->fep_state);
448         }
449
450         if (MASTER(cr))
451         {
452             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
453             if (constructVsites && DOMAINDECOMP(cr))
454             {
455                 gmx_fatal(FARGS,
456                           "Vsite recalculation with -rerun is not implemented with domain "
457                           "decomposition, "
458                           "use a single rank");
459             }
460             if (constructVsites)
461             {
462                 wallcycle_start(wcycle, ewcVSITECONSTR);
463                 vsite->construct(state->x, state->v, state->box, VSiteOperation::PositionsAndVelocities);
464                 wallcycle_stop(wcycle, ewcVSITECONSTR);
465             }
466         }
467
468         if (DOMAINDECOMP(cr))
469         {
470             /* Repartition the domain decomposition */
471             const bool bMasterState = true;
472             dd_partition_system(fplog,
473                                 mdlog,
474                                 step,
475                                 cr,
476                                 bMasterState,
477                                 nstglobalcomm,
478                                 state_global,
479                                 top_global,
480                                 *ir,
481                                 imdSession,
482                                 pull_work,
483                                 state,
484                                 &f,
485                                 mdAtoms,
486                                 &top,
487                                 fr,
488                                 vsite,
489                                 constr,
490                                 nrnb,
491                                 wcycle,
492                                 mdrunOptions.verbose);
493         }
494
495         if (MASTER(cr))
496         {
497             EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
498         }
499
500         if (ir->efep != FreeEnergyPerturbationType::No)
501         {
502             update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
503         }
504
505         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
506                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 is solved
507                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
508
509         if (shellfc)
510         {
511             /* Now is the time to relax the shells */
512             relax_shell_flexcon(fplog,
513                                 cr,
514                                 ms,
515                                 mdrunOptions.verbose,
516                                 enforcedRotation,
517                                 step,
518                                 ir,
519                                 imdSession,
520                                 pull_work,
521                                 bNS,
522                                 force_flags,
523                                 &top,
524                                 constr,
525                                 enerd,
526                                 state->natoms,
527                                 state->x.arrayRefWithPadding(),
528                                 state->v.arrayRefWithPadding(),
529                                 state->box,
530                                 state->lambda,
531                                 &state->hist,
532                                 &f.view(),
533                                 force_vir,
534                                 *mdatoms,
535                                 nrnb,
536                                 wcycle,
537                                 shellfc,
538                                 fr,
539                                 runScheduleWork,
540                                 t,
541                                 mu_tot,
542                                 vsite,
543                                 ddBalanceRegionHandler);
544         }
545         else
546         {
547             /* The coordinates (x) are shifted (to get whole molecules)
548              * in do_force.
549              * This is parallellized as well, and does communication too.
550              * Check comments in sim_util.c
551              */
552             Awh*       awh = nullptr;
553             gmx_edsam* ed  = nullptr;
554             do_force(fplog,
555                      cr,
556                      ms,
557                      *ir,
558                      awh,
559                      enforcedRotation,
560                      imdSession,
561                      pull_work,
562                      step,
563                      nrnb,
564                      wcycle,
565                      &top,
566                      state->box,
567                      state->x.arrayRefWithPadding(),
568                      &state->hist,
569                      &f.view(),
570                      force_vir,
571                      mdatoms,
572                      enerd,
573                      state->lambda,
574                      fr,
575                      runScheduleWork,
576                      vsite,
577                      mu_tot,
578                      t,
579                      ed,
580                      GMX_FORCE_NS | force_flags,
581                      ddBalanceRegionHandler);
582         }
583
584         /* Now we have the energies and forces corresponding to the
585          * coordinates at time t.
586          */
587         {
588             const bool isCheckpointingStep = false;
589             const bool doRerun             = false;
590             const bool bSumEkinhOld        = false;
591             do_md_trajectory_writing(fplog,
592                                      cr,
593                                      nfile,
594                                      fnm,
595                                      step,
596                                      step_rel,
597                                      t,
598                                      ir,
599                                      state,
600                                      state_global,
601                                      observablesHistory,
602                                      top_global,
603                                      fr,
604                                      outf,
605                                      energyOutput,
606                                      ekind,
607                                      f.view().force(),
608                                      isCheckpointingStep,
609                                      doRerun,
610                                      isLastStep,
611                                      mdrunOptions.writeConfout,
612                                      bSumEkinhOld);
613         }
614
615         stopHandler->setSignal();
616
617         {
618             const bool          doInterSimSignal = false;
619             const bool          doIntraSimSignal = true;
620             bool                bSumEkinhOld     = false;
621             t_vcm*              vcm              = nullptr;
622             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
623
624             int cglo_flags = CGLO_GSTAT | CGLO_ENERGY;
625             if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
626             {
627                 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
628             }
629             compute_globals(gstat,
630                             cr,
631                             ir,
632                             fr,
633                             ekind,
634                             makeConstArrayRef(state->x),
635                             makeConstArrayRef(state->v),
636                             state->box,
637                             mdatoms,
638                             nrnb,
639                             vcm,
640                             wcycle,
641                             enerd,
642                             nullptr,
643                             nullptr,
644                             nullptr,
645                             nullptr,
646                             constr != nullptr ? constr->rmsdData() : gmx::ArrayRef<real>{},
647                             &signaller,
648                             state->box,
649                             &bSumEkinhOld,
650                             cglo_flags);
651             if (DOMAINDECOMP(cr))
652             {
653                 checkNumberOfBondedInteractions(
654                         mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
655             }
656         }
657
658         {
659             gmx::HostVector<gmx::RVec>     fglobal(top_global.natoms);
660             gmx::ArrayRef<gmx::RVec>       ftemp;
661             gmx::ArrayRef<const gmx::RVec> flocal = f.view().force();
662             if (DOMAINDECOMP(cr))
663             {
664                 ftemp = gmx::makeArrayRef(fglobal);
665                 dd_collect_vec(cr->dd, state->ddp_count, state->ddp_count_cg_gl, state->cg_gl, flocal, ftemp);
666             }
667             else
668             {
669                 ftemp = f.view().force();
670             }
671
672             if (MASTER(cr))
673             {
674                 MimicCommunicator::sendEnergies(enerd->term[F_EPOT]);
675                 MimicCommunicator::sendForces(ftemp, state_global->natoms);
676             }
677         }
678
679
680         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
681            the virial that should probably be addressed eventually. state->veta has better properies,
682            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
683            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
684
685         if (ir->efep != FreeEnergyPerturbationType::No)
686         {
687             /* Sum up the foreign energy and dhdl terms for md and sd.
688                Currently done every step so that dhdl is correct in the .edr */
689             accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
690         }
691
692         /* Output stuff */
693         if (MASTER(cr))
694         {
695             const bool bCalcEnerStep = true;
696             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
697                                              bCalcEnerStep,
698                                              t,
699                                              mdatoms->tmass,
700                                              enerd,
701                                              ir->fepvals.get(),
702                                              ir->expandedvals.get(),
703                                              state->box,
704                                              PTCouplingArrays({ state->boxv,
705                                                                 state->nosehoover_xi,
706                                                                 state->nosehoover_vxi,
707                                                                 state->nhpres_xi,
708                                                                 state->nhpres_vxi }),
709                                              state->fep_state,
710                                              shake_vir,
711                                              force_vir,
712                                              total_vir,
713                                              pres,
714                                              ekind,
715                                              mu_tot,
716                                              constr);
717
718             const bool do_ene = true;
719             const bool do_log = true;
720             Awh*       awh    = nullptr;
721             const bool do_dr  = ir->nstdisreout != 0;
722             const bool do_or  = ir->nstorireout != 0;
723
724             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
725             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
726                                                do_ene,
727                                                do_dr,
728                                                do_or,
729                                                do_log ? fplog : nullptr,
730                                                step,
731                                                t,
732                                                fr->fcdata.get(),
733                                                awh);
734
735             if (do_per_step(step, ir->nstlog))
736             {
737                 if (fflush(fplog) != 0)
738                 {
739                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
740                 }
741             }
742         }
743
744         /* Print the remaining wall clock time for the run */
745         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
746         {
747             if (shellfc)
748             {
749                 fprintf(stderr, "\n");
750             }
751             print_time(stderr, walltime_accounting, step, ir, cr);
752         }
753
754         cycles = wallcycle_stop(wcycle, ewcSTEP);
755         if (DOMAINDECOMP(cr) && wcycle)
756         {
757             dd_cycles_add(cr->dd, cycles, ddCyclStep);
758         }
759
760         /* increase the MD step number */
761         step++;
762         step_rel++;
763     }
764     /* End of main MD loop */
765
766     /* Closing TNG files can include compressing data. Therefore it is good to do that
767      * before stopping the time measurements. */
768     mdoutf_tng_close(outf);
769
770     /* Stop measuring walltime */
771     walltime_accounting_end_time(walltime_accounting);
772
773     if (MASTER(cr))
774     {
775         MimicCommunicator::finalize();
776     }
777
778     if (!thisRankHasDuty(cr, DUTY_PME))
779     {
780         /* Tell the PME only node to finish */
781         gmx_pme_send_finish(cr);
782     }
783
784     done_mdoutf(outf);
785
786     done_shellfc(fplog, shellfc, step_rel);
787
788     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
789 }