Improve docs and naming for MdModulesNotifiers
[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,
232                        *ir,
233                        gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
234                        MASTER(cr),
235                        &state_global->fep_state,
236                        state_global->lambda);
237
238     const bool        simulationsShareState = false;
239     gmx_mdoutf*       outf                  = init_mdoutf(fplog,
240                                    nfile,
241                                    fnm,
242                                    mdrunOptions,
243                                    cr,
244                                    outputProvider,
245                                    mdModulesNotifiers,
246                                    ir,
247                                    top_global,
248                                    oenv,
249                                    wcycle,
250                                    StartingBehavior::NewSimulation,
251                                    simulationsShareState,
252                                    ms);
253     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
254                                    top_global,
255                                    *ir,
256                                    pull_work,
257                                    mdoutf_get_fp_dhdl(outf),
258                                    true,
259                                    StartingBehavior::NewSimulation,
260                                    simulationsShareState,
261                                    mdModulesNotifiers);
262
263     gstat = global_stat_init(ir);
264
265     /* Check for polarizable models and flexible constraints */
266     shellfc = init_shell_flexcon(fplog,
267                                  top_global,
268                                  constr ? constr->numFlexibleConstraints() : 0,
269                                  ir->nstcalcenergy,
270                                  DOMAINDECOMP(cr),
271                                  runScheduleWork->simulationWork.useGpuPme);
272
273     {
274         double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
275         if ((io > 2000) && MASTER(cr))
276         {
277             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
278         }
279     }
280
281     // Local state only becomes valid now.
282     std::unique_ptr<t_state> stateInstance;
283     t_state*                 state;
284
285     gmx_localtop_t top(top_global.ffparams);
286
287     if (DOMAINDECOMP(cr))
288     {
289         stateInstance = std::make_unique<t_state>();
290         state         = stateInstance.get();
291         dd_init_local_state(*cr->dd, state_global, state);
292
293         /* Distribute the charge groups over the nodes from the master node */
294         dd_partition_system(fplog,
295                             mdlog,
296                             ir->init_step,
297                             cr,
298                             TRUE,
299                             1,
300                             state_global,
301                             top_global,
302                             *ir,
303                             imdSession,
304                             pull_work,
305                             state,
306                             &f,
307                             mdAtoms,
308                             &top,
309                             fr,
310                             vsite,
311                             constr,
312                             nrnb,
313                             nullptr,
314                             FALSE);
315     }
316     else
317     {
318         state_change_natoms(state_global, state_global->natoms);
319         /* Copy the pointer to the global state */
320         state = state_global;
321
322         mdAlgorithmsSetupAtomData(cr, *ir, top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
323     }
324
325     auto mdatoms = mdAtoms->mdatoms();
326
327     // NOTE: The global state is no longer used at this point.
328     // But state_global is still used as temporary storage space for writing
329     // the global state to file and potentially for replica exchange.
330     // (Global topology should persist.)
331
332     update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
333
334     if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
335     {
336         doFreeEnergyPerturbation = true;
337     }
338
339     {
340         int cglo_flags = CGLO_GSTAT;
341         if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
342         {
343             cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
344         }
345         bool   bSumEkinhOld = false;
346         t_vcm* vcm          = nullptr;
347         compute_globals(gstat,
348                         cr,
349                         ir,
350                         fr,
351                         ekind,
352                         makeConstArrayRef(state->x),
353                         makeConstArrayRef(state->v),
354                         state->box,
355                         mdatoms,
356                         nrnb,
357                         vcm,
358                         nullptr,
359                         enerd,
360                         force_vir,
361                         shake_vir,
362                         total_vir,
363                         pres,
364                         gmx::ArrayRef<real>{},
365                         &nullSignaller,
366                         state->box,
367                         &bSumEkinhOld,
368                         cglo_flags);
369         if (DOMAINDECOMP(cr))
370         {
371             checkNumberOfBondedInteractions(
372                     mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
373         }
374     }
375
376     if (MASTER(cr))
377     {
378         fprintf(stderr, "starting MiMiC MD run '%s'\n\n", *(top_global.name));
379         if (mdrunOptions.verbose)
380         {
381             fprintf(stderr,
382                     "Calculated time to finish depends on nsteps from "
383                     "run input file,\nwhich may not correspond to the time "
384                     "needed to process input trajectory.\n\n");
385         }
386         fprintf(fplog, "\n");
387     }
388
389     walltime_accounting_start_time(walltime_accounting);
390     wallcycle_start(wcycle, ewcRUN);
391     print_start(fplog, cr, walltime_accounting, "mdrun");
392
393     /***********************************************************
394      *
395      *             Loop over MD steps
396      *
397      ************************************************************/
398
399     if (constr)
400     {
401         GMX_LOG(mdlog.info)
402                 .asParagraph()
403                 .appendText(
404                         "Simulations has constraints. Constraints will "
405                         "be handled by CPMD.");
406     }
407
408     GMX_LOG(mdlog.info)
409             .asParagraph()
410             .appendText(
411                     "MiMiC does not report kinetic energy, total energy, temperature, virial and "
412                     "pressure.");
413
414     step     = ir->init_step;
415     step_rel = 0;
416
417     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
418             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
419             false,
420             MASTER(cr),
421             ir->nstlist,
422             mdrunOptions.reproducible,
423             nstglobalcomm,
424             mdrunOptions.maximumHoursToRun,
425             ir->nstlist == 0,
426             fplog,
427             step,
428             bNS,
429             walltime_accounting);
430
431     // we don't do counter resetting in rerun - finish will always be valid
432     walltime_accounting_set_valid_finish(walltime_accounting);
433
434     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
435
436     /* and stop now if we should */
437     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
438     while (!isLastStep)
439     {
440         isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel == ir->nsteps));
441         wallcycle_start(wcycle, ewcSTEP);
442
443         t = step;
444
445         if (MASTER(cr))
446         {
447             MimicCommunicator::getCoords(&state_global->x, state_global->natoms);
448         }
449
450         if (ir->efep != FreeEnergyPerturbationType::No)
451         {
452             state->lambda = currentLambdas(step, *(ir->fepvals), state_global->fep_state);
453         }
454
455         if (MASTER(cr))
456         {
457             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
458             if (constructVsites && DOMAINDECOMP(cr))
459             {
460                 gmx_fatal(FARGS,
461                           "Vsite recalculation with -rerun is not implemented with domain "
462                           "decomposition, "
463                           "use a single rank");
464             }
465             if (constructVsites)
466             {
467                 wallcycle_start(wcycle, ewcVSITECONSTR);
468                 vsite->construct(state->x, state->v, state->box, VSiteOperation::PositionsAndVelocities);
469                 wallcycle_stop(wcycle, ewcVSITECONSTR);
470             }
471         }
472
473         if (DOMAINDECOMP(cr))
474         {
475             /* Repartition the domain decomposition */
476             const bool bMasterState = true;
477             dd_partition_system(fplog,
478                                 mdlog,
479                                 step,
480                                 cr,
481                                 bMasterState,
482                                 nstglobalcomm,
483                                 state_global,
484                                 top_global,
485                                 *ir,
486                                 imdSession,
487                                 pull_work,
488                                 state,
489                                 &f,
490                                 mdAtoms,
491                                 &top,
492                                 fr,
493                                 vsite,
494                                 constr,
495                                 nrnb,
496                                 wcycle,
497                                 mdrunOptions.verbose);
498         }
499
500         if (MASTER(cr))
501         {
502             EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
503         }
504
505         if (ir->efep != FreeEnergyPerturbationType::No)
506         {
507             update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
508         }
509
510         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
511                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 is solved
512                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
513
514         if (shellfc)
515         {
516             /* Now is the time to relax the shells */
517             relax_shell_flexcon(fplog,
518                                 cr,
519                                 ms,
520                                 mdrunOptions.verbose,
521                                 enforcedRotation,
522                                 step,
523                                 ir,
524                                 imdSession,
525                                 pull_work,
526                                 bNS,
527                                 force_flags,
528                                 &top,
529                                 constr,
530                                 enerd,
531                                 state->natoms,
532                                 state->x.arrayRefWithPadding(),
533                                 state->v.arrayRefWithPadding(),
534                                 state->box,
535                                 state->lambda,
536                                 &state->hist,
537                                 &f.view(),
538                                 force_vir,
539                                 *mdatoms,
540                                 nrnb,
541                                 wcycle,
542                                 shellfc,
543                                 fr,
544                                 runScheduleWork,
545                                 t,
546                                 mu_tot,
547                                 vsite,
548                                 ddBalanceRegionHandler);
549         }
550         else
551         {
552             /* The coordinates (x) are shifted (to get whole molecules)
553              * in do_force.
554              * This is parallellized as well, and does communication too.
555              * Check comments in sim_util.c
556              */
557             Awh*       awh = nullptr;
558             gmx_edsam* ed  = nullptr;
559             do_force(fplog,
560                      cr,
561                      ms,
562                      *ir,
563                      awh,
564                      enforcedRotation,
565                      imdSession,
566                      pull_work,
567                      step,
568                      nrnb,
569                      wcycle,
570                      &top,
571                      state->box,
572                      state->x.arrayRefWithPadding(),
573                      &state->hist,
574                      &f.view(),
575                      force_vir,
576                      mdatoms,
577                      enerd,
578                      state->lambda,
579                      fr,
580                      runScheduleWork,
581                      vsite,
582                      mu_tot,
583                      t,
584                      ed,
585                      GMX_FORCE_NS | force_flags,
586                      ddBalanceRegionHandler);
587         }
588
589         /* Now we have the energies and forces corresponding to the
590          * coordinates at time t.
591          */
592         {
593             const bool isCheckpointingStep = false;
594             const bool doRerun             = false;
595             const bool bSumEkinhOld        = false;
596             do_md_trajectory_writing(fplog,
597                                      cr,
598                                      nfile,
599                                      fnm,
600                                      step,
601                                      step_rel,
602                                      t,
603                                      ir,
604                                      state,
605                                      state_global,
606                                      observablesHistory,
607                                      top_global,
608                                      fr,
609                                      outf,
610                                      energyOutput,
611                                      ekind,
612                                      f.view().force(),
613                                      isCheckpointingStep,
614                                      doRerun,
615                                      isLastStep,
616                                      mdrunOptions.writeConfout,
617                                      bSumEkinhOld);
618         }
619
620         stopHandler->setSignal();
621
622         {
623             const bool          doInterSimSignal = false;
624             const bool          doIntraSimSignal = true;
625             bool                bSumEkinhOld     = false;
626             t_vcm*              vcm              = nullptr;
627             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
628
629             int cglo_flags = CGLO_GSTAT | CGLO_ENERGY;
630             if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
631             {
632                 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
633             }
634             compute_globals(gstat,
635                             cr,
636                             ir,
637                             fr,
638                             ekind,
639                             makeConstArrayRef(state->x),
640                             makeConstArrayRef(state->v),
641                             state->box,
642                             mdatoms,
643                             nrnb,
644                             vcm,
645                             wcycle,
646                             enerd,
647                             nullptr,
648                             nullptr,
649                             nullptr,
650                             nullptr,
651                             constr != nullptr ? constr->rmsdData() : gmx::ArrayRef<real>{},
652                             &signaller,
653                             state->box,
654                             &bSumEkinhOld,
655                             cglo_flags);
656             if (DOMAINDECOMP(cr))
657             {
658                 checkNumberOfBondedInteractions(
659                         mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
660             }
661         }
662
663         {
664             gmx::HostVector<gmx::RVec>     fglobal(top_global.natoms);
665             gmx::ArrayRef<gmx::RVec>       ftemp;
666             gmx::ArrayRef<const gmx::RVec> flocal = f.view().force();
667             if (DOMAINDECOMP(cr))
668             {
669                 ftemp = gmx::makeArrayRef(fglobal);
670                 dd_collect_vec(cr->dd, state->ddp_count, state->ddp_count_cg_gl, state->cg_gl, flocal, ftemp);
671             }
672             else
673             {
674                 ftemp = f.view().force();
675             }
676
677             if (MASTER(cr))
678             {
679                 MimicCommunicator::sendEnergies(enerd->term[F_EPOT]);
680                 MimicCommunicator::sendForces(ftemp, state_global->natoms);
681             }
682         }
683
684
685         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
686            the virial that should probably be addressed eventually. state->veta has better properies,
687            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
688            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
689
690         if (ir->efep != FreeEnergyPerturbationType::No)
691         {
692             /* Sum up the foreign energy and dhdl terms for md and sd.
693                Currently done every step so that dhdl is correct in the .edr */
694             accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
695         }
696
697         /* Output stuff */
698         if (MASTER(cr))
699         {
700             const bool bCalcEnerStep = true;
701             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
702                                              bCalcEnerStep,
703                                              t,
704                                              mdatoms->tmass,
705                                              enerd,
706                                              ir->fepvals.get(),
707                                              ir->expandedvals.get(),
708                                              state->box,
709                                              PTCouplingArrays({ state->boxv,
710                                                                 state->nosehoover_xi,
711                                                                 state->nosehoover_vxi,
712                                                                 state->nhpres_xi,
713                                                                 state->nhpres_vxi }),
714                                              state->fep_state,
715                                              shake_vir,
716                                              force_vir,
717                                              total_vir,
718                                              pres,
719                                              ekind,
720                                              mu_tot,
721                                              constr);
722
723             const bool do_ene = true;
724             const bool do_log = true;
725             Awh*       awh    = nullptr;
726             const bool do_dr  = ir->nstdisreout != 0;
727             const bool do_or  = ir->nstorireout != 0;
728
729             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
730             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
731                                                do_ene,
732                                                do_dr,
733                                                do_or,
734                                                do_log ? fplog : nullptr,
735                                                step,
736                                                t,
737                                                fr->fcdata.get(),
738                                                awh);
739
740             if (do_per_step(step, ir->nstlog))
741             {
742                 if (fflush(fplog) != 0)
743                 {
744                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
745                 }
746             }
747         }
748
749         /* Print the remaining wall clock time for the run */
750         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
751         {
752             if (shellfc)
753             {
754                 fprintf(stderr, "\n");
755             }
756             print_time(stderr, walltime_accounting, step, ir, cr);
757         }
758
759         cycles = wallcycle_stop(wcycle, ewcSTEP);
760         if (DOMAINDECOMP(cr) && wcycle)
761         {
762             dd_cycles_add(cr->dd, cycles, ddCyclStep);
763         }
764
765         /* increase the MD step number */
766         step++;
767         step_rel++;
768     }
769     /* End of main MD loop */
770
771     /* Closing TNG files can include compressing data. Therefore it is good to do that
772      * before stopping the time measurements. */
773     mdoutf_tng_close(outf);
774
775     /* Stop measuring walltime */
776     walltime_accounting_end_time(walltime_accounting);
777
778     if (MASTER(cr))
779     {
780         MimicCommunicator::finalize();
781     }
782
783     if (!thisRankHasDuty(cr, DUTY_PME))
784     {
785         /* Tell the PME only node to finish */
786         gmx_pme_send_finish(cr);
787     }
788
789     done_mdoutf(outf);
790
791     done_shellfc(fplog, shellfc, step_rel);
792
793     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
794 }