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