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