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