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