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