Refactor md_enums
[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/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/listed_forces.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/freeenergyparameters.h"
85 #include "gromacs/mdlib/md_support.h"
86 #include "gromacs/mdlib/mdatoms.h"
87 #include "gromacs/mdlib/mdoutf.h"
88 #include "gromacs/mdlib/membed.h"
89 #include "gromacs/mdlib/resethandler.h"
90 #include "gromacs/mdlib/sighandler.h"
91 #include "gromacs/mdlib/simulationsignal.h"
92 #include "gromacs/mdlib/stat.h"
93 #include "gromacs/mdlib/stophandler.h"
94 #include "gromacs/mdlib/tgroup.h"
95 #include "gromacs/mdlib/trajectory_writing.h"
96 #include "gromacs/mdlib/update.h"
97 #include "gromacs/mdlib/vcm.h"
98 #include "gromacs/mdlib/vsite.h"
99 #include "gromacs/mdrunutility/handlerestart.h"
100 #include "gromacs/mdrunutility/multisim.h"
101 #include "gromacs/mdrunutility/printtime.h"
102 #include "gromacs/mdtypes/awh_history.h"
103 #include "gromacs/mdtypes/awh_params.h"
104 #include "gromacs/mdtypes/commrec.h"
105 #include "gromacs/mdtypes/df_history.h"
106 #include "gromacs/mdtypes/enerdata.h"
107 #include "gromacs/mdtypes/energyhistory.h"
108 #include "gromacs/mdtypes/forcebuffers.h"
109 #include "gromacs/mdtypes/forcerec.h"
110 #include "gromacs/mdtypes/group.h"
111 #include "gromacs/mdtypes/inputrec.h"
112 #include "gromacs/mdtypes/interaction_const.h"
113 #include "gromacs/mdtypes/md_enums.h"
114 #include "gromacs/mdtypes/mdatom.h"
115 #include "gromacs/mdtypes/mdrunoptions.h"
116 #include "gromacs/mdtypes/observableshistory.h"
117 #include "gromacs/mdtypes/simulation_workload.h"
118 #include "gromacs/mdtypes/state.h"
119 #include "gromacs/mimic/communicator.h"
120 #include "gromacs/mimic/utilities.h"
121 #include "gromacs/pbcutil/pbc.h"
122 #include "gromacs/pulling/pull.h"
123 #include "gromacs/timing/wallcycle.h"
124 #include "gromacs/timing/walltime_accounting.h"
125 #include "gromacs/topology/atoms.h"
126 #include "gromacs/topology/idef.h"
127 #include "gromacs/topology/mtop_util.h"
128 #include "gromacs/topology/topology.h"
129 #include "gromacs/trajectory/trajectoryframe.h"
130 #include "gromacs/utility/basedefinitions.h"
131 #include "gromacs/utility/cstringutil.h"
132 #include "gromacs/utility/fatalerror.h"
133 #include "gromacs/utility/logger.h"
134 #include "gromacs/utility/real.h"
135
136 #include "legacysimulator.h"
137 #include "replicaexchange.h"
138 #include "shellfc.h"
139
140 using gmx::SimulationSignaller;
141
142 void gmx::LegacySimulator::do_mimic()
143 {
144     const t_inputrec* ir = inputrec;
145     int64_t           step, step_rel;
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     /* Domain decomposition could incorrectly miss a bonded
159        interaction, but checking for that requires a global
160        communication stage, which does not otherwise happen in DD
161        code. So we do that alongside the first global energy reduction
162        after a new DD is made. These variables handle whether the
163        check happens, and the result it returns. */
164     bool shouldCheckNumberOfBondedInteractions = false;
165     int  totalNumberOfBondedInteractions       = -1;
166
167     SimulationSignals signals;
168     // Most global communnication stages don't propagate mdrun
169     // signals, and will use this object to achieve that.
170     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
171
172     if (ir->bExpanded)
173     {
174         gmx_fatal(FARGS, "Expanded ensemble not supported by MiMiC.");
175     }
176     if (ir->bSimTemp)
177     {
178         gmx_fatal(FARGS, "Simulated tempering not supported by MiMiC.");
179     }
180     if (ir->bDoAwh)
181     {
182         gmx_fatal(FARGS, "AWH not supported by MiMiC.");
183     }
184     if (replExParams.exchangeInterval > 0)
185     {
186         gmx_fatal(FARGS, "Replica exchange not supported by MiMiC.");
187     }
188     if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
189     {
190         gmx_fatal(FARGS, "Essential dynamics not supported by MiMiC.");
191     }
192     if (ir->bIMD)
193     {
194         gmx_fatal(FARGS, "Interactive MD not supported by MiMiC.");
195     }
196     if (isMultiSim(ms))
197     {
198         gmx_fatal(FARGS, "Multiple simulations not supported by MiMiC.");
199     }
200     if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc, [](SimulatedAnnealing i) {
201             return i != SimulatedAnnealing::No;
202         }))
203     {
204         gmx_fatal(FARGS, "Simulated annealing not supported by MiMiC.");
205     }
206
207     /* Settings for rerun */
208     {
209         // TODO: Avoid changing inputrec (#3854)
210         auto* nonConstInputrec               = const_cast<t_inputrec*>(inputrec);
211         nonConstInputrec->nstlist            = 1;
212         nonConstInputrec->nstcalcenergy      = 1;
213         nonConstInputrec->nstxout_compressed = 0;
214     }
215     int        nstglobalcomm = 1;
216     const bool bNS           = true;
217
218     if (MASTER(cr))
219     {
220         MimicCommunicator::init();
221         auto nonConstGlobalTopology = const_cast<gmx_mtop_t*>(top_global);
222         MimicCommunicator::sendInitData(nonConstGlobalTopology, state_global->x);
223         // TODO: Avoid changing inputrec (#3854)
224         auto* nonConstInputrec   = const_cast<t_inputrec*>(inputrec);
225         nonConstInputrec->nsteps = MimicCommunicator::getStepNumber();
226     }
227     if (DOMAINDECOMP(cr))
228     {
229         // TODO: Avoid changing inputrec (#3854)
230         auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
231         gmx_bcast(sizeof(ir->nsteps), &nonConstInputrec->nsteps, cr->mpi_comm_mygroup);
232     }
233
234     const SimulationGroups* groups = &top_global->groups;
235     {
236         auto nonConstGlobalTopology                          = const_cast<gmx_mtop_t*>(top_global);
237         nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(*top_global);
238     }
239
240     initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda);
241
242     const bool        simulationsShareState = false;
243     gmx_mdoutf*       outf                  = init_mdoutf(fplog,
244                                    nfile,
245                                    fnm,
246                                    mdrunOptions,
247                                    cr,
248                                    outputProvider,
249                                    mdModulesNotifier,
250                                    ir,
251                                    top_global,
252                                    oenv,
253                                    wcycle,
254                                    StartingBehavior::NewSimulation,
255                                    simulationsShareState,
256                                    ms);
257     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
258                                    top_global,
259                                    ir,
260                                    pull_work,
261                                    mdoutf_get_fp_dhdl(outf),
262                                    true,
263                                    StartingBehavior::NewSimulation,
264                                    simulationsShareState,
265                                    mdModulesNotifier);
266
267     gstat = global_stat_init(ir);
268
269     /* Check for polarizable models and flexible constraints */
270     shellfc = init_shell_flexcon(fplog,
271                                  top_global,
272                                  constr ? constr->numFlexibleConstraints() : 0,
273                                  ir->nstcalcenergy,
274                                  DOMAINDECOMP(cr),
275                                  runScheduleWork->simulationWork.useGpuPme);
276
277     {
278         double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
279         if ((io > 2000) && MASTER(cr))
280         {
281             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
282         }
283     }
284
285     // Local state only becomes valid now.
286     std::unique_ptr<t_state> stateInstance;
287     t_state*                 state;
288
289     gmx_localtop_t top(top_global->ffparams);
290
291     if (DOMAINDECOMP(cr))
292     {
293         stateInstance = std::make_unique<t_state>();
294         state         = stateInstance.get();
295         dd_init_local_state(cr->dd, state_global, state);
296
297         /* Distribute the charge groups over the nodes from the master node */
298         dd_partition_system(fplog,
299                             mdlog,
300                             ir->init_step,
301                             cr,
302                             TRUE,
303                             1,
304                             state_global,
305                             *top_global,
306                             ir,
307                             imdSession,
308                             pull_work,
309                             state,
310                             &f,
311                             mdAtoms,
312                             &top,
313                             fr,
314                             vsite,
315                             constr,
316                             nrnb,
317                             nullptr,
318                             FALSE);
319         shouldCheckNumberOfBondedInteractions = true;
320     }
321     else
322     {
323         state_change_natoms(state_global, state_global->natoms);
324         /* Copy the pointer to the global state */
325         state = state_global;
326
327         mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
328     }
329
330     auto mdatoms = mdAtoms->mdatoms();
331
332     // NOTE: The global state is no longer used at this point.
333     // But state_global is still used as temporary storage space for writing
334     // the global state to file and potentially for replica exchange.
335     // (Global topology should persist.)
336
337     update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
338
339     if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
340     {
341         doFreeEnergyPerturbation = true;
342     }
343
344     {
345         int cglo_flags =
346                 (CGLO_GSTAT
347                  | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
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                         &totalNumberOfBondedInteractions,
371                         &bSumEkinhOld,
372                         cglo_flags);
373     }
374     checkNumberOfBondedInteractions(mdlog,
375                                     cr,
376                                     totalNumberOfBondedInteractions,
377                                     top_global,
378                                     &top,
379                                     makeConstArrayRef(state->x),
380                                     state->box,
381                                     &shouldCheckNumberOfBondedInteractions);
382
383     if (MASTER(cr))
384     {
385         fprintf(stderr, "starting MiMiC MD run '%s'\n\n", *(top_global->name));
386         if (mdrunOptions.verbose)
387         {
388             fprintf(stderr,
389                     "Calculated time to finish depends on nsteps from "
390                     "run input file,\nwhich may not correspond to the time "
391                     "needed to process input trajectory.\n\n");
392         }
393         fprintf(fplog, "\n");
394     }
395
396     walltime_accounting_start_time(walltime_accounting);
397     wallcycle_start(wcycle, ewcRUN);
398     print_start(fplog, cr, walltime_accounting, "mdrun");
399
400     /***********************************************************
401      *
402      *             Loop over MD steps
403      *
404      ************************************************************/
405
406     if (constr)
407     {
408         GMX_LOG(mdlog.info)
409                 .asParagraph()
410                 .appendText(
411                         "Simulations has constraints. Constraints will "
412                         "be handled by CPMD.");
413     }
414
415     GMX_LOG(mdlog.info)
416             .asParagraph()
417             .appendText(
418                     "MiMiC does not report kinetic energy, total energy, temperature, virial and "
419                     "pressure.");
420
421     step     = ir->init_step;
422     step_rel = 0;
423
424     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
425             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
426             false,
427             MASTER(cr),
428             ir->nstlist,
429             mdrunOptions.reproducible,
430             nstglobalcomm,
431             mdrunOptions.maximumHoursToRun,
432             ir->nstlist == 0,
433             fplog,
434             step,
435             bNS,
436             walltime_accounting);
437
438     // we don't do counter resetting in rerun - finish will always be valid
439     walltime_accounting_set_valid_finish(walltime_accounting);
440
441     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
442
443     /* and stop now if we should */
444     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
445     while (!isLastStep)
446     {
447         isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel == ir->nsteps));
448         wallcycle_start(wcycle, ewcSTEP);
449
450         t = step;
451
452         if (MASTER(cr))
453         {
454             MimicCommunicator::getCoords(&state_global->x, state_global->natoms);
455         }
456
457         if (ir->efep != FreeEnergyPerturbationType::No)
458         {
459             state->lambda = currentLambdas(step, *(ir->fepvals), state_global->fep_state);
460         }
461
462         if (MASTER(cr))
463         {
464             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
465             if (constructVsites && DOMAINDECOMP(cr))
466             {
467                 gmx_fatal(FARGS,
468                           "Vsite recalculation with -rerun is not implemented with domain "
469                           "decomposition, "
470                           "use a single rank");
471             }
472             if (constructVsites)
473             {
474                 wallcycle_start(wcycle, ewcVSITECONSTR);
475                 vsite->construct(state->x, state->v, state->box, VSiteOperation::PositionsAndVelocities);
476                 wallcycle_stop(wcycle, ewcVSITECONSTR);
477             }
478         }
479
480         if (DOMAINDECOMP(cr))
481         {
482             /* Repartition the domain decomposition */
483             const bool bMasterState = true;
484             dd_partition_system(fplog,
485                                 mdlog,
486                                 step,
487                                 cr,
488                                 bMasterState,
489                                 nstglobalcomm,
490                                 state_global,
491                                 *top_global,
492                                 ir,
493                                 imdSession,
494                                 pull_work,
495                                 state,
496                                 &f,
497                                 mdAtoms,
498                                 &top,
499                                 fr,
500                                 vsite,
501                                 constr,
502                                 nrnb,
503                                 wcycle,
504                                 mdrunOptions.verbose);
505             shouldCheckNumberOfBondedInteractions = true;
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             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                             &totalNumberOfBondedInteractions,
658                             &bSumEkinhOld,
659                             CGLO_GSTAT | CGLO_ENERGY
660                                     | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
661                                                                              : 0));
662             checkNumberOfBondedInteractions(mdlog,
663                                             cr,
664                                             totalNumberOfBondedInteractions,
665                                             top_global,
666                                             &top,
667                                             makeConstArrayRef(state->x),
668                                             state->box,
669                                             &shouldCheckNumberOfBondedInteractions);
670         }
671
672         {
673             gmx::HostVector<gmx::RVec>     fglobal(top_global->natoms);
674             gmx::ArrayRef<gmx::RVec>       ftemp;
675             gmx::ArrayRef<const gmx::RVec> flocal = f.view().force();
676             if (DOMAINDECOMP(cr))
677             {
678                 ftemp = gmx::makeArrayRef(fglobal);
679                 dd_collect_vec(cr->dd, state->ddp_count, state->ddp_count_cg_gl, state->cg_gl, flocal, ftemp);
680             }
681             else
682             {
683                 ftemp = f.view().force();
684             }
685
686             if (MASTER(cr))
687             {
688                 MimicCommunicator::sendEnergies(enerd->term[F_EPOT]);
689                 MimicCommunicator::sendForces(ftemp, state_global->natoms);
690             }
691         }
692
693
694         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
695            the virial that should probably be addressed eventually. state->veta has better properies,
696            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
697            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
698
699         if (ir->efep != FreeEnergyPerturbationType::No)
700         {
701             /* Sum up the foreign energy and dhdl terms for md and sd.
702                Currently done every step so that dhdl is correct in the .edr */
703             accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
704         }
705
706         /* Output stuff */
707         if (MASTER(cr))
708         {
709             const bool bCalcEnerStep = true;
710             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
711                                              bCalcEnerStep,
712                                              t,
713                                              mdatoms->tmass,
714                                              enerd,
715                                              ir->fepvals.get(),
716                                              ir->expandedvals.get(),
717                                              state->box,
718                                              PTCouplingArrays({ state->boxv,
719                                                                 state->nosehoover_xi,
720                                                                 state->nosehoover_vxi,
721                                                                 state->nhpres_xi,
722                                                                 state->nhpres_vxi }),
723                                              state->fep_state,
724                                              shake_vir,
725                                              force_vir,
726                                              total_vir,
727                                              pres,
728                                              ekind,
729                                              mu_tot,
730                                              constr);
731
732             const bool do_ene = true;
733             const bool do_log = true;
734             Awh*       awh    = nullptr;
735             const bool do_dr  = ir->nstdisreout != 0;
736             const bool do_or  = ir->nstorireout != 0;
737
738             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
739             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
740                                                do_ene,
741                                                do_dr,
742                                                do_or,
743                                                do_log ? fplog : nullptr,
744                                                step,
745                                                t,
746                                                fr->fcdata.get(),
747                                                awh);
748
749             if (do_per_step(step, ir->nstlog))
750             {
751                 if (fflush(fplog) != 0)
752                 {
753                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
754                 }
755             }
756         }
757
758         /* Print the remaining wall clock time for the run */
759         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
760         {
761             if (shellfc)
762             {
763                 fprintf(stderr, "\n");
764             }
765             print_time(stderr, walltime_accounting, step, ir, cr);
766         }
767
768         cycles = wallcycle_stop(wcycle, ewcSTEP);
769         if (DOMAINDECOMP(cr) && wcycle)
770         {
771             dd_cycles_add(cr->dd, cycles, ddCyclStep);
772         }
773
774         /* increase the MD step number */
775         step++;
776         step_rel++;
777     }
778     /* End of main MD loop */
779
780     /* Closing TNG files can include compressing data. Therefore it is good to do that
781      * before stopping the time measurements. */
782     mdoutf_tng_close(outf);
783
784     /* Stop measuring walltime */
785     walltime_accounting_end_time(walltime_accounting);
786
787     if (MASTER(cr))
788     {
789         MimicCommunicator::finalize();
790     }
791
792     if (!thisRankHasDuty(cr, DUTY_PME))
793     {
794         /* Tell the PME only node to finish */
795         gmx_pme_send_finish(cr);
796     }
797
798     done_mdoutf(outf);
799
800     done_shellfc(fplog, shellfc, step_rel);
801
802     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
803 }