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