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