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