Merge "Merge release-2019 into master"
[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, 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/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.h"
64 #include "gromacs/ewald/pme_load_balancing.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/manage_threading.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/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/mdoutf.h"
87 #include "gromacs/mdlib/membed.h"
88 #include "gromacs/mdlib/resethandler.h"
89 #include "gromacs/mdlib/sighandler.h"
90 #include "gromacs/mdlib/simulationsignal.h"
91 #include "gromacs/mdlib/stat.h"
92 #include "gromacs/mdlib/stophandler.h"
93 #include "gromacs/mdlib/tgroup.h"
94 #include "gromacs/mdlib/trajectory_writing.h"
95 #include "gromacs/mdlib/update.h"
96 #include "gromacs/mdlib/vcm.h"
97 #include "gromacs/mdlib/vsite.h"
98 #include "gromacs/mdrunutility/handlerestart.h"
99 #include "gromacs/mdrunutility/multisim.h"
100 #include "gromacs/mdrunutility/printtime.h"
101 #include "gromacs/mdtypes/awh_history.h"
102 #include "gromacs/mdtypes/awh_params.h"
103 #include "gromacs/mdtypes/commrec.h"
104 #include "gromacs/mdtypes/df_history.h"
105 #include "gromacs/mdtypes/enerdata.h"
106 #include "gromacs/mdtypes/energyhistory.h"
107 #include "gromacs/mdtypes/fcdata.h"
108 #include "gromacs/mdtypes/forcerec.h"
109 #include "gromacs/mdtypes/group.h"
110 #include "gromacs/mdtypes/inputrec.h"
111 #include "gromacs/mdtypes/interaction_const.h"
112 #include "gromacs/mdtypes/md_enums.h"
113 #include "gromacs/mdtypes/mdatom.h"
114 #include "gromacs/mdtypes/mdrunoptions.h"
115 #include "gromacs/mdtypes/observableshistory.h"
116 #include "gromacs/mdtypes/state.h"
117 #include "gromacs/mimic/communicator.h"
118 #include "gromacs/mimic/utilities.h"
119 #include "gromacs/pbcutil/mshift.h"
120 #include "gromacs/pbcutil/pbc.h"
121 #include "gromacs/pulling/pull.h"
122 #include "gromacs/timing/wallcycle.h"
123 #include "gromacs/timing/walltime_accounting.h"
124 #include "gromacs/topology/atoms.h"
125 #include "gromacs/topology/idef.h"
126 #include "gromacs/topology/mtop_util.h"
127 #include "gromacs/topology/topology.h"
128 #include "gromacs/trajectory/trajectoryframe.h"
129 #include "gromacs/utility/basedefinitions.h"
130 #include "gromacs/utility/cstringutil.h"
131 #include "gromacs/utility/fatalerror.h"
132 #include "gromacs/utility/logger.h"
133 #include "gromacs/utility/real.h"
134 #include "gromacs/utility/smalloc.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, lam0[efptNR];
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     gmx_localtop_t           top;
153     PaddedVector<gmx::RVec>  f {};
154     gmx_global_stat_t        gstat;
155     t_graph                 *graph = nullptr;
156     gmx_shellfc_t           *shellfc;
157
158     double                   cycles;
159
160     /* Domain decomposition could incorrectly miss a bonded
161        interaction, but checking for that requires a global
162        communication stage, which does not otherwise happen in DD
163        code. So we do that alongside the first global energy reduction
164        after a new DD is made. These variables handle whether the
165        check happens, and the result it returns. */
166     bool              shouldCheckNumberOfBondedInteractions = false;
167     int               totalNumberOfBondedInteractions       = -1;
168
169     SimulationSignals signals;
170     // Most global communnication stages don't propagate mdrun
171     // signals, and will use this object to achieve that.
172     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
173
174     if (ir->bExpanded)
175     {
176         gmx_fatal(FARGS, "Expanded ensemble not supported by MiMiC.");
177     }
178     if (ir->bSimTemp)
179     {
180         gmx_fatal(FARGS, "Simulated tempering not supported by MiMiC.");
181     }
182     if (ir->bDoAwh)
183     {
184         gmx_fatal(FARGS, "AWH not supported by MiMiC.");
185     }
186     if (replExParams.exchangeInterval > 0)
187     {
188         gmx_fatal(FARGS, "Replica exchange not supported by MiMiC.");
189     }
190     if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
191     {
192         gmx_fatal(FARGS, "Essential dynamics not supported by MiMiC.");
193     }
194     if (ir->bIMD)
195     {
196         gmx_fatal(FARGS, "Interactive MD not supported by MiMiC.");
197     }
198     if (isMultiSim(ms))
199     {
200         gmx_fatal(FARGS, "Multiple simulations not supported by MiMiC.");
201     }
202     if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
203                     [](int i){return i != eannNO; }))
204     {
205         gmx_fatal(FARGS, "Simulated annealing not supported by MiMiC.");
206     }
207
208     /* Settings for rerun */
209     ir->nstlist       = 1;
210     ir->nstcalcenergy = 1;
211     int        nstglobalcomm = 1;
212     const bool bNS           = true;
213
214     // Communicator to interact with MiMiC
215     MimicCommunicator mimicCommunicator {};
216     if (MASTER(cr))
217     {
218         mimicCommunicator.init();
219         mimicCommunicator.sendInitData(top_global, state_global->x);
220         ir->nsteps = mimicCommunicator.getStepNumber();
221     }
222
223     ir->nstxout_compressed                   = 0;
224     SimulationGroups *groups                        = &top_global->groups;
225     top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
226
227     initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
228
229     gmx_mdoutf        *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, oenv, wcycle,
230                                           StartingBehavior::NewSimulation);
231     gmx::EnergyOutput  energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work, mdoutf_get_fp_dhdl(outf), true);
232
233     /* Kinetic energy data */
234     std::unique_ptr<gmx_ekindata_t> eKinData = std::make_unique<gmx_ekindata_t>();
235     gmx_ekindata_t                 *ekind    = eKinData.get();
236     init_ekindata(fplog, top_global, &(ir->opts), ekind);
237     /* Copy the cos acceleration to the groups struct */
238     ekind->cosacc.cos_accel = ir->cos_accel;
239
240     gstat = global_stat_init(ir);
241
242     /* Check for polarizable models and flexible constraints */
243     shellfc = init_shell_flexcon(fplog,
244                                  top_global, constr ? constr->numFlexibleConstraints() : 0,
245                                  ir->nstcalcenergy, DOMAINDECOMP(cr));
246
247     {
248         double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
249         if ((io > 2000) && MASTER(cr))
250         {
251             fprintf(stderr,
252                     "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
253                     io);
254         }
255     }
256
257     // Local state only becomes valid now.
258     std::unique_ptr<t_state> stateInstance;
259     t_state *                state;
260
261     if (DOMAINDECOMP(cr))
262     {
263         dd_init_local_top(*top_global, &top);
264
265         stateInstance = std::make_unique<t_state>();
266         state         = stateInstance.get();
267         dd_init_local_state(cr->dd, state_global, state);
268
269         /* Distribute the charge groups over the nodes from the master node */
270         dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
271                             state_global, *top_global, ir, imdSession,
272                             pull_work,
273                             state, &f, mdAtoms, &top, fr,
274                             vsite, constr,
275                             nrnb, nullptr, FALSE);
276         shouldCheckNumberOfBondedInteractions = true;
277         gmx_bcast(sizeof(ir->nsteps), &ir->nsteps, cr);
278     }
279     else
280     {
281         state_change_natoms(state_global, state_global->natoms);
282         /* We need to allocate one element extra, since we might use
283          * (unaligned) 4-wide SIMD loads to access rvec entries.
284          */
285         f.resizeWithPadding(state_global->natoms);
286         /* Copy the pointer to the global state */
287         state = state_global;
288
289         mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr,
290                                   &graph, mdAtoms, constr, vsite, shellfc);
291     }
292
293     auto mdatoms = mdAtoms->mdatoms();
294
295     // NOTE: The global state is no longer used at this point.
296     // But state_global is still used as temporary storage space for writing
297     // the global state to file and potentially for replica exchange.
298     // (Global topology should persist.)
299
300     update_mdatoms(mdatoms, state->lambda[efptMASS]);
301
302     if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
303     {
304         doFreeEnergyPerturbation = true;
305     }
306
307     {
308         int    cglo_flags = (CGLO_GSTAT |
309                              (shouldCheckNumberOfBondedInteractions ?
310                               CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
311         bool   bSumEkinhOld = false;
312         t_vcm *vcm          = nullptr;
313         compute_globals(gstat, cr, ir, fr, ekind,
314                         state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
315                         mdatoms, nrnb, vcm,
316                         nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
317                         constr, &nullSignaller, state->box,
318                         &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
319     }
320     checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
321                                     top_global, &top, state->x.rvec_array(), state->box,
322                                     &shouldCheckNumberOfBondedInteractions);
323
324     if (MASTER(cr))
325     {
326         fprintf(stderr, "starting MiMiC MD run '%s'\n\n",
327                 *(top_global->name));
328         if (mdrunOptions.verbose)
329         {
330             fprintf(stderr, "Calculated time to finish depends on nsteps from "
331                     "run input file,\nwhich may not correspond to the time "
332                     "needed to process input trajectory.\n\n");
333         }
334         fprintf(fplog, "\n");
335     }
336
337     walltime_accounting_start_time(walltime_accounting);
338     wallcycle_start(wcycle, ewcRUN);
339     print_start(fplog, cr, walltime_accounting, "mdrun");
340
341     /***********************************************************
342      *
343      *             Loop over MD steps
344      *
345      ************************************************************/
346
347     if (constr)
348     {
349         GMX_LOG(mdlog.info).asParagraph().
350             appendText("Simulations has constraints. Constraints will "
351                        "be handled by CPMD.");
352     }
353
354     GMX_LOG(mdlog.info).asParagraph().
355         appendText("MiMiC does not report kinetic energy, total energy, temperature, virial and pressure.");
356
357     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
358                 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false,
359                 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstglobalcomm,
360                 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
361
362     // we don't do counter resetting in rerun - finish will always be valid
363     walltime_accounting_set_valid_finish(walltime_accounting);
364
365     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
366
367     step     = ir->init_step;
368     step_rel = 0;
369
370     /* and stop now if we should */
371     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
372     while (!isLastStep)
373     {
374         isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel == ir->nsteps));
375         wallcycle_start(wcycle, ewcSTEP);
376
377         t = step;
378
379         if (MASTER(cr))
380         {
381             mimicCommunicator.getCoords(&state_global->x, state_global->natoms);
382         }
383
384         if (ir->efep != efepNO)
385         {
386             setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
387         }
388
389         if (MASTER(cr))
390         {
391             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
392             if (constructVsites && DOMAINDECOMP(cr))
393             {
394                 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, "
395                           "use a single rank");
396             }
397         }
398
399         if (DOMAINDECOMP(cr))
400         {
401             /* Repartition the domain decomposition */
402             const bool bMasterState = true;
403             dd_partition_system(fplog, mdlog, step, cr,
404                                 bMasterState, nstglobalcomm,
405                                 state_global, *top_global, ir, imdSession,
406                                 pull_work,
407                                 state, &f, mdAtoms, &top, fr,
408                                 vsite, constr,
409                                 nrnb, wcycle,
410                                 mdrunOptions.verbose);
411             shouldCheckNumberOfBondedInteractions = true;
412         }
413
414         if (MASTER(cr))
415         {
416             energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
417         }
418
419         if (ir->efep != efepNO)
420         {
421             update_mdatoms(mdatoms, state->lambda[efptMASS]);
422         }
423
424         force_flags = (GMX_FORCE_STATECHANGED |
425                        GMX_FORCE_DYNAMICBOX |
426                        GMX_FORCE_ALLFORCES |
427                        GMX_FORCE_VIRIAL |  // TODO: Get rid of this once #2649 is solved
428                        GMX_FORCE_ENERGY |
429                        (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
430
431         if (shellfc)
432         {
433             /* Now is the time to relax the shells */
434             relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
435                                 enforcedRotation, step,
436                                 ir, imdSession, pull_work, bNS, force_flags, &top,
437                                 constr, enerd, fcd,
438                                 state->natoms,
439                                 state->x.arrayRefWithPadding(),
440                                 state->v.arrayRefWithPadding(),
441                                 state->box,
442                                 state->lambda,
443                                 &state->hist,
444                                 f.arrayRefWithPadding(), force_vir, mdatoms,
445                                 nrnb, wcycle, graph,
446                                 shellfc, fr, ppForceWorkload, t, mu_tot,
447                                 vsite,
448                                 ddBalanceRegionHandler);
449         }
450         else
451         {
452             /* The coordinates (x) are shifted (to get whole molecules)
453              * in do_force.
454              * This is parallellized as well, and does communication too.
455              * Check comments in sim_util.c
456              */
457             Awh       *awh = nullptr;
458             gmx_edsam *ed  = nullptr;
459             do_force(fplog, cr, ms, ir, awh, enforcedRotation, imdSession,
460                      pull_work,
461                      step, nrnb, wcycle, &top,
462                      state->box, state->x.arrayRefWithPadding(), &state->hist,
463                      f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
464                      state->lambda, graph,
465                      fr, ppForceWorkload, vsite, mu_tot, t, ed,
466                      GMX_FORCE_NS | force_flags,
467                      ddBalanceRegionHandler);
468         }
469
470         /* Now we have the energies and forces corresponding to the
471          * coordinates at time t.
472          */
473         {
474             const bool isCheckpointingStep = false;
475             const bool doRerun             = false;
476             const bool bSumEkinhOld        = false;
477             do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
478                                      ir, state, state_global, observablesHistory,
479                                      top_global, fr,
480                                      outf, energyOutput, ekind, f,
481                                      isCheckpointingStep, doRerun, isLastStep,
482                                      mdrunOptions.writeConfout,
483                                      bSumEkinhOld);
484         }
485
486         stopHandler->setSignal();
487
488         if (graph)
489         {
490             /* Need to unshift here */
491             unshift_self(graph, state->box, as_rvec_array(state->x.data()));
492         }
493
494         if (vsite != nullptr)
495         {
496             wallcycle_start(wcycle, ewcVSITECONSTR);
497             if (graph != nullptr)
498             {
499                 shift_self(graph, state->box, as_rvec_array(state->x.data()));
500             }
501             construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
502                              top.idef.iparams, top.idef.il,
503                              fr->ePBC, fr->bMolPBC, cr, state->box);
504
505             if (graph != nullptr)
506             {
507                 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
508             }
509             wallcycle_stop(wcycle, ewcVSITECONSTR);
510         }
511
512         {
513             const bool          doInterSimSignal = false;
514             const bool          doIntraSimSignal = true;
515             bool                bSumEkinhOld     = false;
516             t_vcm              *vcm              = nullptr;
517             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
518
519             compute_globals(gstat, cr, ir, fr, ekind,
520                             state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
521                             mdatoms, nrnb, vcm,
522                             wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
523                             constr, &signaller,
524                             state->box,
525                             &totalNumberOfBondedInteractions, &bSumEkinhOld,
526                             CGLO_GSTAT | CGLO_ENERGY
527                             | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
528                             );
529             checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
530                                             top_global, &top, state->x.rvec_array(), state->box,
531                                             &shouldCheckNumberOfBondedInteractions);
532         }
533
534         {
535             gmx::HostVector<gmx::RVec>     fglobal(top_global->natoms);
536             gmx::ArrayRef<gmx::RVec>       ftemp;
537             gmx::ArrayRef<const gmx::RVec> flocal = gmx::makeArrayRef(f);
538             if (DOMAINDECOMP(cr))
539             {
540                 ftemp = gmx::makeArrayRef(fglobal);
541                 dd_collect_vec(cr->dd, state, flocal, ftemp);
542             }
543             else
544             {
545                 ftemp = gmx::makeArrayRef(f);
546             }
547
548             if (MASTER(cr))
549             {
550                 mimicCommunicator.sendEnergies(enerd->term[F_EPOT]);
551                 mimicCommunicator.sendForces(ftemp, state_global->natoms);
552             }
553         }
554
555
556
557         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
558            the virial that should probably be addressed eventually. state->veta has better properies,
559            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
560            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
561
562         if (ir->efep != efepNO)
563         {
564             /* Sum up the foreign energy and dhdl terms for md and sd.
565                Currently done every step so that dhdl is correct in the .edr */
566             sum_dhdl(enerd, state->lambda, ir->fepvals);
567         }
568
569         /* Output stuff */
570         if (MASTER(cr))
571         {
572             const bool bCalcEnerStep = true;
573             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep,
574                                              t, mdatoms->tmass, enerd, state,
575                                              ir->fepvals, ir->expandedvals, state->box,
576                                              shake_vir, force_vir, total_vir, pres,
577                                              ekind, mu_tot, constr);
578
579             const bool do_ene = true;
580             const bool do_log = true;
581             Awh       *awh    = nullptr;
582             const bool do_dr  = ir->nstdisreout != 0;
583             const bool do_or  = ir->nstorireout != 0;
584
585             energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
586             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
587                                                do_log ? fplog : nullptr,
588                                                step, t,
589                                                fcd, awh);
590
591             if (do_per_step(step, ir->nstlog))
592             {
593                 if (fflush(fplog) != 0)
594                 {
595                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
596                 }
597             }
598         }
599
600         /* Print the remaining wall clock time for the run */
601         if (isMasterSimMasterRank(ms, MASTER(cr)) &&
602             (mdrunOptions.verbose || gmx_got_usr_signal()))
603         {
604             if (shellfc)
605             {
606                 fprintf(stderr, "\n");
607             }
608             print_time(stderr, walltime_accounting, step, ir, cr);
609         }
610
611         cycles = wallcycle_stop(wcycle, ewcSTEP);
612         if (DOMAINDECOMP(cr) && wcycle)
613         {
614             dd_cycles_add(cr->dd, cycles, ddCyclStep);
615         }
616
617         /* increase the MD step number */
618         step++;
619         step_rel++;
620     }
621     /* End of main MD loop */
622
623     /* Closing TNG files can include compressing data. Therefore it is good to do that
624      * before stopping the time measurements. */
625     mdoutf_tng_close(outf);
626
627     /* Stop measuring walltime */
628     walltime_accounting_end_time(walltime_accounting);
629
630     if (MASTER(cr))
631     {
632         mimicCommunicator.finalize();
633     }
634
635     if (!thisRankHasDuty(cr, DUTY_PME))
636     {
637         /* Tell the PME only node to finish */
638         gmx_pme_send_finish(cr);
639     }
640
641     done_mdoutf(outf);
642
643     done_shellfc(fplog, shellfc, step_rel);
644
645     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
646 }