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