Restrict scope of step and step_rel
[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,2021, 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/applied_forces/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/localtopologychecker.h"
61 #include "gromacs/domdec/mdsetup.h"
62 #include "gromacs/domdec/partition.h"
63 #include "gromacs/essentialdynamics/edsam.h"
64 #include "gromacs/ewald/pme_load_balancing.h"
65 #include "gromacs/ewald/pme_pp.h"
66 #include "gromacs/fileio/trxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/gpu_utils/gpu_utils.h"
70 #include "gromacs/listed_forces/listed_forces.h"
71 #include "gromacs/math/functions.h"
72 #include "gromacs/math/utilities.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/math/vectypes.h"
75 #include "gromacs/mdlib/checkpointhandler.h"
76 #include "gromacs/mdlib/compute_io.h"
77 #include "gromacs/mdlib/constr.h"
78 #include "gromacs/mdlib/ebin.h"
79 #include "gromacs/mdlib/enerdata_utils.h"
80 #include "gromacs/mdlib/energyoutput.h"
81 #include "gromacs/mdlib/expanded.h"
82 #include "gromacs/mdlib/force.h"
83 #include "gromacs/mdlib/force_flags.h"
84 #include "gromacs/mdlib/forcerec.h"
85 #include "gromacs/mdlib/freeenergyparameters.h"
86 #include "gromacs/mdlib/md_support.h"
87 #include "gromacs/mdlib/mdatoms.h"
88 #include "gromacs/mdlib/mdoutf.h"
89 #include "gromacs/mdlib/membed.h"
90 #include "gromacs/mdlib/resethandler.h"
91 #include "gromacs/mdlib/sighandler.h"
92 #include "gromacs/mdlib/simulationsignal.h"
93 #include "gromacs/mdlib/stat.h"
94 #include "gromacs/mdlib/stophandler.h"
95 #include "gromacs/mdlib/tgroup.h"
96 #include "gromacs/mdlib/trajectory_writing.h"
97 #include "gromacs/mdlib/update.h"
98 #include "gromacs/mdlib/vcm.h"
99 #include "gromacs/mdlib/vsite.h"
100 #include "gromacs/mdrunutility/handlerestart.h"
101 #include "gromacs/mdrunutility/multisim.h"
102 #include "gromacs/mdrunutility/printtime.h"
103 #include "gromacs/mdtypes/awh_history.h"
104 #include "gromacs/mdtypes/awh_params.h"
105 #include "gromacs/mdtypes/commrec.h"
106 #include "gromacs/mdtypes/df_history.h"
107 #include "gromacs/mdtypes/enerdata.h"
108 #include "gromacs/mdtypes/energyhistory.h"
109 #include "gromacs/mdtypes/forcebuffers.h"
110 #include "gromacs/mdtypes/forcerec.h"
111 #include "gromacs/mdtypes/group.h"
112 #include "gromacs/mdtypes/inputrec.h"
113 #include "gromacs/mdtypes/interaction_const.h"
114 #include "gromacs/mdtypes/md_enums.h"
115 #include "gromacs/mdtypes/mdatom.h"
116 #include "gromacs/mdtypes/mdrunoptions.h"
117 #include "gromacs/mdtypes/observableshistory.h"
118 #include "gromacs/mdtypes/simulation_workload.h"
119 #include "gromacs/mdtypes/state.h"
120 #include "gromacs/mimic/communicator.h"
121 #include "gromacs/mimic/utilities.h"
122 #include "gromacs/pbcutil/pbc.h"
123 #include "gromacs/pulling/pull.h"
124 #include "gromacs/timing/wallcycle.h"
125 #include "gromacs/timing/walltime_accounting.h"
126 #include "gromacs/topology/atoms.h"
127 #include "gromacs/topology/idef.h"
128 #include "gromacs/topology/mtop_util.h"
129 #include "gromacs/topology/topology.h"
130 #include "gromacs/trajectory/trajectoryframe.h"
131 #include "gromacs/utility/basedefinitions.h"
132 #include "gromacs/utility/cstringutil.h"
133 #include "gromacs/utility/fatalerror.h"
134 #include "gromacs/utility/logger.h"
135 #include "gromacs/utility/real.h"
136
137 #include "legacysimulator.h"
138 #include "replicaexchange.h"
139 #include "shellfc.h"
140
141 using gmx::SimulationSignaller;
142
143 void gmx::LegacySimulator::do_mimic()
144 {
145     const t_inputrec* ir = inputrec;
146     double            t;
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     ForceBuffers      f;
153     gmx_global_stat_t gstat;
154     gmx_shellfc_t*    shellfc;
155
156     double cycles;
157
158     SimulationSignals signals;
159     // Most global communnication stages don't propagate mdrun
160     // signals, and will use this object to achieve that.
161     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
162
163     if (ir->bExpanded)
164     {
165         gmx_fatal(FARGS, "Expanded ensemble not supported by MiMiC.");
166     }
167     if (ir->bSimTemp)
168     {
169         gmx_fatal(FARGS, "Simulated tempering not supported by MiMiC.");
170     }
171     if (ir->bDoAwh)
172     {
173         gmx_fatal(FARGS, "AWH not supported by MiMiC.");
174     }
175     if (replExParams.exchangeInterval > 0)
176     {
177         gmx_fatal(FARGS, "Replica exchange not supported by MiMiC.");
178     }
179     if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
180     {
181         gmx_fatal(FARGS, "Essential dynamics not supported by MiMiC.");
182     }
183     if (ir->bIMD)
184     {
185         gmx_fatal(FARGS, "Interactive MD not supported by MiMiC.");
186     }
187     if (isMultiSim(ms))
188     {
189         gmx_fatal(FARGS, "Multiple simulations not supported by MiMiC.");
190     }
191     if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc, [](SimulatedAnnealing i) {
192             return i != SimulatedAnnealing::No;
193         }))
194     {
195         gmx_fatal(FARGS, "Simulated annealing not supported by MiMiC.");
196     }
197
198     /* Settings for rerun */
199     {
200         // TODO: Avoid changing inputrec (#3854)
201         auto* nonConstInputrec               = const_cast<t_inputrec*>(inputrec);
202         nonConstInputrec->nstlist            = 1;
203         nonConstInputrec->nstcalcenergy      = 1;
204         nonConstInputrec->nstxout_compressed = 0;
205     }
206     int        nstglobalcomm = 1;
207     const bool bNS           = true;
208
209     if (MASTER(cr))
210     {
211         MimicCommunicator::init();
212         auto* nonConstGlobalTopology = const_cast<gmx_mtop_t*>(&top_global);
213         MimicCommunicator::sendInitData(nonConstGlobalTopology, state_global->x);
214         // TODO: Avoid changing inputrec (#3854)
215         auto* nonConstInputrec   = const_cast<t_inputrec*>(inputrec);
216         nonConstInputrec->nsteps = MimicCommunicator::getStepNumber();
217     }
218     if (DOMAINDECOMP(cr))
219     {
220         // TODO: Avoid changing inputrec (#3854)
221         auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
222         gmx_bcast(sizeof(ir->nsteps), &nonConstInputrec->nsteps, cr->mpi_comm_mygroup);
223     }
224
225     const SimulationGroups* groups = &top_global.groups;
226     {
227         auto* nonConstGlobalTopology                         = const_cast<gmx_mtop_t*>(&top_global);
228         nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(top_global);
229     }
230
231     initialize_lambdas(fplog,
232                        ir->efep,
233                        ir->bSimTemp,
234                        *ir->fepvals,
235                        ir->simtempvals->temperatures,
236                        gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
237                        MASTER(cr),
238                        &state_global->fep_state,
239                        state_global->lambda);
240
241     const bool        simulationsShareState = false;
242     gmx_mdoutf*       outf                  = init_mdoutf(fplog,
243                                    nfile,
244                                    fnm,
245                                    mdrunOptions,
246                                    cr,
247                                    outputProvider,
248                                    mdModulesNotifiers,
249                                    ir,
250                                    top_global,
251                                    oenv,
252                                    wcycle,
253                                    StartingBehavior::NewSimulation,
254                                    simulationsShareState,
255                                    ms);
256     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
257                                    top_global,
258                                    *ir,
259                                    pull_work,
260                                    mdoutf_get_fp_dhdl(outf),
261                                    true,
262                                    StartingBehavior::NewSimulation,
263                                    simulationsShareState,
264                                    mdModulesNotifiers);
265
266     gstat = global_stat_init(ir);
267
268     /* Check for polarizable models and flexible constraints */
269     shellfc = init_shell_flexcon(fplog,
270                                  top_global,
271                                  constr ? constr->numFlexibleConstraints() : 0,
272                                  ir->nstcalcenergy,
273                                  DOMAINDECOMP(cr),
274                                  runScheduleWork->simulationWork.useGpuPme);
275
276     {
277         double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
278         if ((io > 2000) && MASTER(cr))
279         {
280             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
281         }
282     }
283
284     // Local state only becomes valid now.
285     std::unique_ptr<t_state> stateInstance;
286     t_state*                 state;
287
288     gmx_localtop_t top(top_global.ffparams);
289
290     if (DOMAINDECOMP(cr))
291     {
292         stateInstance = std::make_unique<t_state>();
293         state         = stateInstance.get();
294         dd_init_local_state(*cr->dd, state_global, state);
295
296         /* Distribute the charge groups over the nodes from the master node */
297         dd_partition_system(fplog,
298                             mdlog,
299                             ir->init_step,
300                             cr,
301                             TRUE,
302                             1,
303                             state_global,
304                             top_global,
305                             *ir,
306                             imdSession,
307                             pull_work,
308                             state,
309                             &f,
310                             mdAtoms,
311                             &top,
312                             fr,
313                             vsite,
314                             constr,
315                             nrnb,
316                             nullptr,
317                             FALSE);
318     }
319     else
320     {
321         state_change_natoms(state_global, state_global->natoms);
322         /* Copy the pointer to the global state */
323         state = state_global;
324
325         mdAlgorithmsSetupAtomData(cr, *ir, top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
326     }
327
328     auto* mdatoms = mdAtoms->mdatoms();
329
330     // NOTE: The global state is no longer used at this point.
331     // But state_global is still used as temporary storage space for writing
332     // the global state to file and potentially for replica exchange.
333     // (Global topology should persist.)
334
335     update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
336
337     if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
338     {
339         doFreeEnergyPerturbation = true;
340     }
341
342     {
343         int cglo_flags = CGLO_GSTAT;
344         if (DOMAINDECOMP(cr) && dd_localTopologyChecker(*cr->dd).shouldCheckNumberOfBondedInteractions())
345         {
346             cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
347         }
348         bool   bSumEkinhOld = false;
349         t_vcm* vcm          = nullptr;
350         compute_globals(gstat,
351                         cr,
352                         ir,
353                         fr,
354                         ekind,
355                         makeConstArrayRef(state->x),
356                         makeConstArrayRef(state->v),
357                         state->box,
358                         mdatoms,
359                         nrnb,
360                         vcm,
361                         nullptr,
362                         enerd,
363                         force_vir,
364                         shake_vir,
365                         total_vir,
366                         pres,
367                         gmx::ArrayRef<real>{},
368                         &nullSignaller,
369                         state->box,
370                         &bSumEkinhOld,
371                         cglo_flags);
372         if (DOMAINDECOMP(cr))
373         {
374             dd_localTopologyChecker(cr->dd)->checkNumberOfBondedInteractions(
375                     &top, makeConstArrayRef(state->x), state->box);
376         }
377     }
378
379     if (MASTER(cr))
380     {
381         fprintf(stderr, "starting MiMiC MD run '%s'\n\n", *(top_global.name));
382         if (mdrunOptions.verbose)
383         {
384             fprintf(stderr,
385                     "Calculated time to finish depends on nsteps from "
386                     "run input file,\nwhich may not correspond to the time "
387                     "needed to process input trajectory.\n\n");
388         }
389         fprintf(fplog, "\n");
390     }
391
392     walltime_accounting_start_time(walltime_accounting);
393     wallcycle_start(wcycle, WallCycleCounter::Run);
394     print_start(fplog, cr, walltime_accounting, "mdrun");
395
396     /***********************************************************
397      *
398      *             Loop over MD steps
399      *
400      ************************************************************/
401
402     if (constr)
403     {
404         GMX_LOG(mdlog.info)
405                 .asParagraph()
406                 .appendText(
407                         "Simulations has constraints. Constraints will "
408                         "be handled by CPMD.");
409     }
410
411     GMX_LOG(mdlog.info)
412             .asParagraph()
413             .appendText(
414                     "MiMiC does not report kinetic energy, total energy, temperature, virial and "
415                     "pressure.");
416
417     int64_t step     = ir->init_step;
418     int64_t step_rel = 0;
419
420     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
421             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
422             false,
423             MASTER(cr),
424             ir->nstlist,
425             mdrunOptions.reproducible,
426             nstglobalcomm,
427             mdrunOptions.maximumHoursToRun,
428             ir->nstlist == 0,
429             fplog,
430             step,
431             bNS,
432             walltime_accounting);
433
434     // we don't do counter resetting in rerun - finish will always be valid
435     walltime_accounting_set_valid_finish(walltime_accounting);
436
437     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
438
439     /* and stop now if we should */
440     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
441     while (!isLastStep)
442     {
443         isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel == ir->nsteps));
444         wallcycle_start(wcycle, WallCycleCounter::Step);
445
446         t = step;
447
448         if (MASTER(cr))
449         {
450             MimicCommunicator::getCoords(&state_global->x, state_global->natoms);
451         }
452
453         if (ir->efep != FreeEnergyPerturbationType::No)
454         {
455             state->lambda = currentLambdas(step, *(ir->fepvals), state_global->fep_state);
456         }
457
458         if (MASTER(cr))
459         {
460             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
461             if (constructVsites && DOMAINDECOMP(cr))
462             {
463                 gmx_fatal(FARGS,
464                           "Vsite recalculation with -rerun is not implemented with domain "
465                           "decomposition, "
466                           "use a single rank");
467             }
468             if (constructVsites)
469             {
470                 wallcycle_start(wcycle, WallCycleCounter::VsiteConstr);
471                 vsite->construct(state->x, state->v, state->box, VSiteOperation::PositionsAndVelocities);
472                 wallcycle_stop(wcycle, WallCycleCounter::VsiteConstr);
473             }
474         }
475
476         if (DOMAINDECOMP(cr))
477         {
478             /* Repartition the domain decomposition */
479             const bool bMasterState = true;
480             dd_partition_system(fplog,
481                                 mdlog,
482                                 step,
483                                 cr,
484                                 bMasterState,
485                                 nstglobalcomm,
486                                 state_global,
487                                 top_global,
488                                 *ir,
489                                 imdSession,
490                                 pull_work,
491                                 state,
492                                 &f,
493                                 mdAtoms,
494                                 &top,
495                                 fr,
496                                 vsite,
497                                 constr,
498                                 nrnb,
499                                 wcycle,
500                                 mdrunOptions.verbose);
501         }
502
503         if (MASTER(cr))
504         {
505             EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
506         }
507
508         if (ir->efep != FreeEnergyPerturbationType::No)
509         {
510             update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
511         }
512
513         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
514                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 is solved
515                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
516
517         if (shellfc)
518         {
519             /* Now is the time to relax the shells */
520             relax_shell_flexcon(fplog,
521                                 cr,
522                                 ms,
523                                 mdrunOptions.verbose,
524                                 enforcedRotation,
525                                 step,
526                                 ir,
527                                 imdSession,
528                                 pull_work,
529                                 bNS,
530                                 force_flags,
531                                 &top,
532                                 constr,
533                                 enerd,
534                                 state->natoms,
535                                 state->x.arrayRefWithPadding(),
536                                 state->v.arrayRefWithPadding(),
537                                 state->box,
538                                 state->lambda,
539                                 &state->hist,
540                                 &f.view(),
541                                 force_vir,
542                                 *mdatoms,
543                                 nrnb,
544                                 wcycle,
545                                 shellfc,
546                                 fr,
547                                 runScheduleWork,
548                                 t,
549                                 mu_tot,
550                                 vsite,
551                                 ddBalanceRegionHandler);
552         }
553         else
554         {
555             /* The coordinates (x) are shifted (to get whole molecules)
556              * in do_force.
557              * This is parallellized as well, and does communication too.
558              * Check comments in sim_util.c
559              */
560             Awh*       awh = nullptr;
561             gmx_edsam* ed  = nullptr;
562             do_force(fplog,
563                      cr,
564                      ms,
565                      *ir,
566                      awh,
567                      enforcedRotation,
568                      imdSession,
569                      pull_work,
570                      step,
571                      nrnb,
572                      wcycle,
573                      &top,
574                      state->box,
575                      state->x.arrayRefWithPadding(),
576                      &state->hist,
577                      &f.view(),
578                      force_vir,
579                      mdatoms,
580                      enerd,
581                      state->lambda,
582                      fr,
583                      runScheduleWork,
584                      vsite,
585                      mu_tot,
586                      t,
587                      ed,
588                      GMX_FORCE_NS | force_flags,
589                      ddBalanceRegionHandler);
590         }
591
592         /* Now we have the energies and forces corresponding to the
593          * coordinates at time t.
594          */
595         {
596             const bool isCheckpointingStep = false;
597             const bool doRerun             = false;
598             const bool bSumEkinhOld        = false;
599             do_md_trajectory_writing(fplog,
600                                      cr,
601                                      nfile,
602                                      fnm,
603                                      step,
604                                      step_rel,
605                                      t,
606                                      ir,
607                                      state,
608                                      state_global,
609                                      observablesHistory,
610                                      top_global,
611                                      fr,
612                                      outf,
613                                      energyOutput,
614                                      ekind,
615                                      f.view().force(),
616                                      isCheckpointingStep,
617                                      doRerun,
618                                      isLastStep,
619                                      mdrunOptions.writeConfout,
620                                      bSumEkinhOld);
621         }
622
623         stopHandler->setSignal();
624
625         {
626             const bool          doInterSimSignal = false;
627             const bool          doIntraSimSignal = true;
628             bool                bSumEkinhOld     = false;
629             t_vcm*              vcm              = nullptr;
630             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
631
632             int cglo_flags = CGLO_GSTAT | CGLO_ENERGY;
633             if (DOMAINDECOMP(cr) && dd_localTopologyChecker(*cr->dd).shouldCheckNumberOfBondedInteractions())
634             {
635                 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
636             }
637             compute_globals(gstat,
638                             cr,
639                             ir,
640                             fr,
641                             ekind,
642                             makeConstArrayRef(state->x),
643                             makeConstArrayRef(state->v),
644                             state->box,
645                             mdatoms,
646                             nrnb,
647                             vcm,
648                             wcycle,
649                             enerd,
650                             nullptr,
651                             nullptr,
652                             nullptr,
653                             nullptr,
654                             constr != nullptr ? constr->rmsdData() : gmx::ArrayRef<real>{},
655                             &signaller,
656                             state->box,
657                             &bSumEkinhOld,
658                             cglo_flags);
659             if (DOMAINDECOMP(cr))
660             {
661                 dd_localTopologyChecker(cr->dd)->checkNumberOfBondedInteractions(
662                         &top, makeConstArrayRef(state->x), state->box);
663             }
664         }
665
666         {
667             gmx::HostVector<gmx::RVec>     fglobal(top_global.natoms);
668             gmx::ArrayRef<gmx::RVec>       ftemp;
669             gmx::ArrayRef<const gmx::RVec> flocal = f.view().force();
670             if (DOMAINDECOMP(cr))
671             {
672                 ftemp = gmx::makeArrayRef(fglobal);
673                 dd_collect_vec(cr->dd, state->ddp_count, state->ddp_count_cg_gl, state->cg_gl, flocal, ftemp);
674             }
675             else
676             {
677                 ftemp = f.view().force();
678             }
679
680             if (MASTER(cr))
681             {
682                 MimicCommunicator::sendEnergies(enerd->term[F_EPOT]);
683                 MimicCommunicator::sendForces(ftemp, state_global->natoms);
684             }
685         }
686
687
688         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
689            the virial that should probably be addressed eventually. state->veta has better properies,
690            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
691            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
692
693         if (ir->efep != FreeEnergyPerturbationType::No)
694         {
695             /* Sum up the foreign energy and dhdl terms for md and sd.
696                Currently done every step so that dhdl is correct in the .edr */
697             accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
698         }
699
700         /* Output stuff */
701         if (MASTER(cr))
702         {
703             const bool bCalcEnerStep = true;
704             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
705                                              bCalcEnerStep,
706                                              t,
707                                              mdatoms->tmass,
708                                              enerd,
709                                              ir->fepvals.get(),
710                                              ir->expandedvals.get(),
711                                              state->box,
712                                              PTCouplingArrays({ state->boxv,
713                                                                 state->nosehoover_xi,
714                                                                 state->nosehoover_vxi,
715                                                                 state->nhpres_xi,
716                                                                 state->nhpres_vxi }),
717                                              state->fep_state,
718                                              total_vir,
719                                              pres,
720                                              ekind,
721                                              mu_tot,
722                                              constr);
723
724             const bool do_ene = true;
725             const bool do_log = true;
726             Awh*       awh    = nullptr;
727             const bool do_dr  = ir->nstdisreout != 0;
728             const bool do_or  = ir->nstorireout != 0;
729
730             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
731             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
732                                                do_ene,
733                                                do_dr,
734                                                do_or,
735                                                do_log ? fplog : nullptr,
736                                                step,
737                                                t,
738                                                fr->fcdata.get(),
739                                                awh);
740
741             if (do_per_step(step, ir->nstlog))
742             {
743                 if (fflush(fplog) != 0)
744                 {
745                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
746                 }
747             }
748         }
749
750         /* Print the remaining wall clock time for the run */
751         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
752         {
753             if (shellfc)
754             {
755                 fprintf(stderr, "\n");
756             }
757             print_time(stderr, walltime_accounting, step, ir, cr);
758         }
759
760         cycles = wallcycle_stop(wcycle, WallCycleCounter::Step);
761         if (DOMAINDECOMP(cr) && wcycle)
762         {
763             dd_cycles_add(cr->dd, cycles, ddCyclStep);
764         }
765
766         /* increase the MD step number */
767         step++;
768         step_rel++;
769     }
770     /* End of main MD loop */
771
772     /* Closing TNG files can include compressing data. Therefore it is good to do that
773      * before stopping the time measurements. */
774     mdoutf_tng_close(outf);
775
776     /* Stop measuring walltime */
777     walltime_accounting_end_time(walltime_accounting);
778
779     if (MASTER(cr))
780     {
781         MimicCommunicator::finalize();
782     }
783
784     if (!thisRankHasDuty(cr, DUTY_PME))
785     {
786         /* Tell the PME only node to finish */
787         gmx_pme_send_finish(cr);
788     }
789
790     done_mdoutf(outf);
791
792     done_shellfc(fplog, shellfc, step_rel);
793
794     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
795 }