Change the behavior of the GPU update UI
[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     PaddedHostVector<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, mdModulesNotifier,
230                                    ir, top_global, oenv, wcycle, StartingBehavior::NewSimulation);
231     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
232                                    mdoutf_get_fp_dhdl(outf), true, StartingBehavior::NewSimulation,
233                                    mdModulesNotifier);
234
235     gstat = global_stat_init(ir);
236
237     /* Check for polarizable models and flexible constraints */
238     shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
239                                  ir->nstcalcenergy, DOMAINDECOMP(cr));
240
241     {
242         double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
243         if ((io > 2000) && MASTER(cr))
244         {
245             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
246         }
247     }
248
249     // Local state only becomes valid now.
250     std::unique_ptr<t_state> stateInstance;
251     t_state*                 state;
252
253     if (DOMAINDECOMP(cr))
254     {
255         dd_init_local_top(*top_global, &top);
256
257         stateInstance = std::make_unique<t_state>();
258         state         = stateInstance.get();
259         dd_init_local_state(cr->dd, state_global, state);
260
261         /* Distribute the charge groups over the nodes from the master node */
262         dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
263                             imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
264                             nrnb, nullptr, FALSE);
265         shouldCheckNumberOfBondedInteractions = true;
266         gmx_bcast(sizeof(ir->nsteps), &ir->nsteps, cr);
267     }
268     else
269     {
270         state_change_natoms(state_global, state_global->natoms);
271         /* We need to allocate one element extra, since we might use
272          * (unaligned) 4-wide SIMD loads to access rvec entries.
273          */
274         f.resizeWithPadding(state_global->natoms);
275         /* Copy the pointer to the global state */
276         state = state_global;
277
278         mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &graph, mdAtoms, constr, vsite, shellfc);
279     }
280
281     auto mdatoms = mdAtoms->mdatoms();
282
283     // NOTE: The global state is no longer used at this point.
284     // But state_global is still used as temporary storage space for writing
285     // the global state to file and potentially for replica exchange.
286     // (Global topology should persist.)
287
288     update_mdatoms(mdatoms, state->lambda[efptMASS]);
289
290     if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
291     {
292         doFreeEnergyPerturbation = true;
293     }
294
295     {
296         int cglo_flags =
297                 (CGLO_GSTAT
298                  | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
299         bool   bSumEkinhOld = false;
300         t_vcm* vcm          = nullptr;
301         compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
302                         state->box, state->lambda[efptVDW], mdatoms, nrnb, vcm, nullptr, enerd,
303                         force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
304                         state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
305     }
306     checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
307                                     state->x.rvec_array(), state->box,
308                                     &shouldCheckNumberOfBondedInteractions);
309
310     if (MASTER(cr))
311     {
312         fprintf(stderr, "starting MiMiC MD run '%s'\n\n", *(top_global->name));
313         if (mdrunOptions.verbose)
314         {
315             fprintf(stderr,
316                     "Calculated time to finish depends on nsteps from "
317                     "run input file,\nwhich may not correspond to the time "
318                     "needed to process input trajectory.\n\n");
319         }
320         fprintf(fplog, "\n");
321     }
322
323     walltime_accounting_start_time(walltime_accounting);
324     wallcycle_start(wcycle, ewcRUN);
325     print_start(fplog, cr, walltime_accounting, "mdrun");
326
327     /***********************************************************
328      *
329      *             Loop over MD steps
330      *
331      ************************************************************/
332
333     if (constr)
334     {
335         GMX_LOG(mdlog.info)
336                 .asParagraph()
337                 .appendText(
338                         "Simulations has constraints. Constraints will "
339                         "be handled by CPMD.");
340     }
341
342     GMX_LOG(mdlog.info)
343             .asParagraph()
344             .appendText(
345                     "MiMiC does not report kinetic energy, total energy, temperature, virial and "
346                     "pressure.");
347
348     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
349             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false, MASTER(cr),
350             ir->nstlist, mdrunOptions.reproducible, nstglobalcomm, mdrunOptions.maximumHoursToRun,
351             ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
352
353     // we don't do counter resetting in rerun - finish will always be valid
354     walltime_accounting_set_valid_finish(walltime_accounting);
355
356     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
357
358     step     = ir->init_step;
359     step_rel = 0;
360
361     /* and stop now if we should */
362     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
363     while (!isLastStep)
364     {
365         isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel == ir->nsteps));
366         wallcycle_start(wcycle, ewcSTEP);
367
368         t = step;
369
370         if (MASTER(cr))
371         {
372             mimicCommunicator.getCoords(&state_global->x, state_global->natoms);
373         }
374
375         if (ir->efep != efepNO)
376         {
377             setCurrentLambdasLocal(step, ir->fepvals, lam0, state->lambda, state->fep_state);
378         }
379
380         if (MASTER(cr))
381         {
382             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
383             if (constructVsites && DOMAINDECOMP(cr))
384             {
385                 gmx_fatal(FARGS,
386                           "Vsite recalculation with -rerun is not implemented with domain "
387                           "decomposition, "
388                           "use a single rank");
389             }
390         }
391
392         if (DOMAINDECOMP(cr))
393         {
394             /* Repartition the domain decomposition */
395             const bool bMasterState = true;
396             dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
397                                 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
398                                 fr, vsite, constr, nrnb, wcycle, mdrunOptions.verbose);
399             shouldCheckNumberOfBondedInteractions = true;
400         }
401
402         if (MASTER(cr))
403         {
404             energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
405         }
406
407         if (ir->efep != efepNO)
408         {
409             update_mdatoms(mdatoms, state->lambda[efptMASS]);
410         }
411
412         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
413                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 is solved
414                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
415
416         if (shellfc)
417         {
418             /* Now is the time to relax the shells */
419             relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
420                                 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
421                                 state->natoms, state->x.arrayRefWithPadding(),
422                                 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
423                                 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, graph,
424                                 shellfc, fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
425         }
426         else
427         {
428             /* The coordinates (x) are shifted (to get whole molecules)
429              * in do_force.
430              * This is parallellized as well, and does communication too.
431              * Check comments in sim_util.c
432              */
433             Awh*       awh = nullptr;
434             gmx_edsam* ed  = nullptr;
435             do_force(fplog, cr, ms, ir, awh, enforcedRotation, imdSession, pull_work, step, nrnb,
436                      wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
437                      f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, graph,
438                      fr, runScheduleWork, vsite, mu_tot, t, ed, GMX_FORCE_NS | force_flags,
439                      ddBalanceRegionHandler);
440         }
441
442         /* Now we have the energies and forces corresponding to the
443          * coordinates at time t.
444          */
445         {
446             const bool isCheckpointingStep = false;
447             const bool doRerun             = false;
448             const bool bSumEkinhOld        = false;
449             do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state,
450                                      state_global, observablesHistory, top_global, fr, outf,
451                                      energyOutput, ekind, f, isCheckpointingStep, doRerun,
452                                      isLastStep, mdrunOptions.writeConfout, bSumEkinhOld);
453         }
454
455         stopHandler->setSignal();
456
457         if (graph)
458         {
459             /* Need to unshift here */
460             unshift_self(graph, state->box, as_rvec_array(state->x.data()));
461         }
462
463         if (vsite != nullptr)
464         {
465             wallcycle_start(wcycle, ewcVSITECONSTR);
466             if (graph != nullptr)
467             {
468                 shift_self(graph, state->box, as_rvec_array(state->x.data()));
469             }
470             construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t,
471                              as_rvec_array(state->v.data()), top.idef.iparams, top.idef.il,
472                              fr->ePBC, fr->bMolPBC, cr, state->box);
473
474             if (graph != nullptr)
475             {
476                 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
477             }
478             wallcycle_stop(wcycle, ewcVSITECONSTR);
479         }
480
481         {
482             const bool          doInterSimSignal = false;
483             const bool          doIntraSimSignal = true;
484             bool                bSumEkinhOld     = false;
485             t_vcm*              vcm              = nullptr;
486             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
487
488             compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
489                             state->box, state->lambda[efptVDW], mdatoms, nrnb, vcm, wcycle, enerd,
490                             nullptr, nullptr, nullptr, nullptr, mu_tot, constr, &signaller,
491                             state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
492                             CGLO_GSTAT | CGLO_ENERGY
493                                     | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
494                                                                              : 0));
495             checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
496                                             &top, state->x.rvec_array(), state->box,
497                                             &shouldCheckNumberOfBondedInteractions);
498         }
499
500         {
501             gmx::HostVector<gmx::RVec>     fglobal(top_global->natoms);
502             gmx::ArrayRef<gmx::RVec>       ftemp;
503             gmx::ArrayRef<const gmx::RVec> flocal = gmx::makeArrayRef(f);
504             if (DOMAINDECOMP(cr))
505             {
506                 ftemp = gmx::makeArrayRef(fglobal);
507                 dd_collect_vec(cr->dd, state, flocal, ftemp);
508             }
509             else
510             {
511                 ftemp = gmx::makeArrayRef(f);
512             }
513
514             if (MASTER(cr))
515             {
516                 mimicCommunicator.sendEnergies(enerd->term[F_EPOT]);
517                 mimicCommunicator.sendForces(ftemp, state_global->natoms);
518             }
519         }
520
521
522         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
523            the virial that should probably be addressed eventually. state->veta has better properies,
524            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
525            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
526
527         if (ir->efep != efepNO)
528         {
529             /* Sum up the foreign energy and dhdl terms for md and sd.
530                Currently done every step so that dhdl is correct in the .edr */
531             sum_dhdl(enerd, state->lambda, *ir->fepvals);
532         }
533
534         /* Output stuff */
535         if (MASTER(cr))
536         {
537             const bool bCalcEnerStep = true;
538             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep, t,
539                                              mdatoms->tmass, enerd, state, ir->fepvals,
540                                              ir->expandedvals, state->box, shake_vir, force_vir,
541                                              total_vir, pres, ekind, mu_tot, constr);
542
543             const bool do_ene = true;
544             const bool do_log = true;
545             Awh*       awh    = nullptr;
546             const bool do_dr  = ir->nstdisreout != 0;
547             const bool do_or  = ir->nstorireout != 0;
548
549             energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
550             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
551                                                do_log ? fplog : nullptr, step, t, fcd, awh);
552
553             if (do_per_step(step, ir->nstlog))
554             {
555                 if (fflush(fplog) != 0)
556                 {
557                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
558                 }
559             }
560         }
561
562         /* Print the remaining wall clock time for the run */
563         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
564         {
565             if (shellfc)
566             {
567                 fprintf(stderr, "\n");
568             }
569             print_time(stderr, walltime_accounting, step, ir, cr);
570         }
571
572         cycles = wallcycle_stop(wcycle, ewcSTEP);
573         if (DOMAINDECOMP(cr) && wcycle)
574         {
575             dd_cycles_add(cr->dd, cycles, ddCyclStep);
576         }
577
578         /* increase the MD step number */
579         step++;
580         step_rel++;
581     }
582     /* End of main MD loop */
583
584     /* Closing TNG files can include compressing data. Therefore it is good to do that
585      * before stopping the time measurements. */
586     mdoutf_tng_close(outf);
587
588     /* Stop measuring walltime */
589     walltime_accounting_end_time(walltime_accounting);
590
591     if (MASTER(cr))
592     {
593         mimicCommunicator.finalize();
594     }
595
596     if (!thisRankHasDuty(cr, DUTY_PME))
597     {
598         /* Tell the PME only node to finish */
599         gmx_pme_send_finish(cr);
600     }
601
602     done_mdoutf(outf);
603
604     done_shellfc(fplog, shellfc, step_rel);
605
606     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
607 }