Fix out of sync checkpoint files in simulations sharing state
[alexxy/gromacs.git] / src / gromacs / mdrun / md.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2011-2019,2020, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /*! \internal \file
38  *
39  * \brief Implements the integrator for normal molecular dynamics simulations
40  *
41  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42  * \ingroup module_mdrun
43  */
44 #include "gmxpre.h"
45
46 #include <cinttypes>
47 #include <cmath>
48 #include <cstdio>
49 #include <cstdlib>
50
51 #include <algorithm>
52 #include <memory>
53
54 #include "gromacs/awh/awh.h"
55 #include "gromacs/commandline/filenm.h"
56 #include "gromacs/domdec/collect.h"
57 #include "gromacs/domdec/dlbtiming.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_network.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/domdec/mdsetup.h"
62 #include "gromacs/domdec/partition.h"
63 #include "gromacs/essentialdynamics/edsam.h"
64 #include "gromacs/ewald/pme.h"
65 #include "gromacs/ewald/pme_load_balancing.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/imd/imd.h"
71 #include "gromacs/listed_forces/manage_threading.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/utilities.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/math/vectypes.h"
76 #include "gromacs/mdlib/checkpointhandler.h"
77 #include "gromacs/mdlib/compute_io.h"
78 #include "gromacs/mdlib/constr.h"
79 #include "gromacs/mdlib/ebin.h"
80 #include "gromacs/mdlib/enerdata_utils.h"
81 #include "gromacs/mdlib/energyoutput.h"
82 #include "gromacs/mdlib/expanded.h"
83 #include "gromacs/mdlib/force.h"
84 #include "gromacs/mdlib/force_flags.h"
85 #include "gromacs/mdlib/forcerec.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/update_constrain_cuda.h"
99 #include "gromacs/mdlib/vcm.h"
100 #include "gromacs/mdlib/vsite.h"
101 #include "gromacs/mdrunutility/handlerestart.h"
102 #include "gromacs/mdrunutility/multisim.h"
103 #include "gromacs/mdrunutility/printtime.h"
104 #include "gromacs/mdtypes/awh_history.h"
105 #include "gromacs/mdtypes/awh_params.h"
106 #include "gromacs/mdtypes/commrec.h"
107 #include "gromacs/mdtypes/df_history.h"
108 #include "gromacs/mdtypes/energyhistory.h"
109 #include "gromacs/mdtypes/fcdata.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/pullhistory.h"
119 #include "gromacs/mdtypes/simulation_workload.h"
120 #include "gromacs/mdtypes/state.h"
121 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
122 #include "gromacs/modularsimulator/energyelement.h"
123 #include "gromacs/nbnxm/gpu_data_mgmt.h"
124 #include "gromacs/nbnxm/nbnxm.h"
125 #include "gromacs/pbcutil/mshift.h"
126 #include "gromacs/pbcutil/pbc.h"
127 #include "gromacs/pulling/output.h"
128 #include "gromacs/pulling/pull.h"
129 #include "gromacs/swap/swapcoords.h"
130 #include "gromacs/timing/wallcycle.h"
131 #include "gromacs/timing/walltime_accounting.h"
132 #include "gromacs/topology/atoms.h"
133 #include "gromacs/topology/idef.h"
134 #include "gromacs/topology/mtop_util.h"
135 #include "gromacs/topology/topology.h"
136 #include "gromacs/trajectory/trajectoryframe.h"
137 #include "gromacs/utility/basedefinitions.h"
138 #include "gromacs/utility/cstringutil.h"
139 #include "gromacs/utility/fatalerror.h"
140 #include "gromacs/utility/logger.h"
141 #include "gromacs/utility/real.h"
142 #include "gromacs/utility/smalloc.h"
143
144 #include "legacysimulator.h"
145 #include "replicaexchange.h"
146 #include "shellfc.h"
147
148 #if GMX_FAHCORE
149 #    include "corewrap.h"
150 #endif
151
152 using gmx::SimulationSignaller;
153
154 void gmx::LegacySimulator::do_md()
155 {
156     // TODO Historically, the EM and MD "integrators" used different
157     // names for the t_inputrec *parameter, but these must have the
158     // same name, now that it's a member of a struct. We use this ir
159     // alias to avoid a large ripple of nearly useless changes.
160     // t_inputrec is being replaced by IMdpOptionsProvider, so this
161     // will go away eventually.
162     t_inputrec*  ir = inputrec;
163     int64_t      step, step_rel;
164     double       t, t0 = ir->init_t, lam0[efptNR];
165     gmx_bool     bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
166     gmx_bool     bNS, bNStList, bStopCM, bFirstStep, bInitStep, bLastStep = FALSE;
167     gmx_bool     bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
168     gmx_bool     do_ene, do_log, do_verbose;
169     gmx_bool     bMasterState;
170     unsigned int force_flags;
171     tensor force_vir = { { 0 } }, shake_vir = { { 0 } }, total_vir = { { 0 } }, tmp_vir = { { 0 } },
172            pres = { { 0 } };
173     int                         i, m;
174     rvec                        mu_tot;
175     matrix                      pressureCouplingMu, M;
176     gmx_repl_ex_t               repl_ex = nullptr;
177     gmx_localtop_t              top;
178     PaddedHostVector<gmx::RVec> f{};
179     gmx_global_stat_t           gstat;
180     t_graph*                    graph = nullptr;
181     gmx_shellfc_t*              shellfc;
182     gmx_bool                    bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
183     gmx_bool                    bTemp, bPres, bTrotter;
184     real                        dvdl_constr;
185     std::vector<RVec>           cbuf;
186     matrix                      lastbox;
187     int                         lamnew = 0;
188     /* for FEP */
189     int       nstfep = 0;
190     double    cycles;
191     real      saved_conserved_quantity = 0;
192     real      last_ekin                = 0;
193     t_extmass MassQ;
194     char      sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
195
196     /* PME load balancing data for GPU kernels */
197     gmx_bool bPMETune         = FALSE;
198     gmx_bool bPMETunePrinting = FALSE;
199
200     bool bInteractiveMDstep = false;
201
202     /* Domain decomposition could incorrectly miss a bonded
203        interaction, but checking for that requires a global
204        communication stage, which does not otherwise happen in DD
205        code. So we do that alongside the first global energy reduction
206        after a new DD is made. These variables handle whether the
207        check happens, and the result it returns. */
208     bool shouldCheckNumberOfBondedInteractions = false;
209     int  totalNumberOfBondedInteractions       = -1;
210
211     SimulationSignals signals;
212     // Most global communnication stages don't propagate mdrun
213     // signals, and will use this object to achieve that.
214     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
215
216     if (!mdrunOptions.writeConfout)
217     {
218         // This is on by default, and the main known use case for
219         // turning it off is for convenience in benchmarking, which is
220         // something that should not show up in the general user
221         // interface.
222         GMX_LOG(mdlog.info)
223                 .asParagraph()
224                 .appendText(
225                         "The -noconfout functionality is deprecated, and may be removed in a "
226                         "future version.");
227     }
228
229     /* md-vv uses averaged full step velocities for T-control
230        md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
231        md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
232     bTrotter = (EI_VV(ir->eI)
233                 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
234
235     const bool bRerunMD = false;
236
237     int nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
238     bGStatEveryStep   = (nstglobalcomm == 1);
239
240     SimulationGroups* groups = &top_global->groups;
241
242     std::unique_ptr<EssentialDynamics> ed = nullptr;
243     if (opt2bSet("-ei", nfile, fnm))
244     {
245         /* Initialize essential dynamics sampling */
246         ed = init_edsam(mdlog, opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm), top_global,
247                         ir, cr, constr, state_global, observablesHistory, oenv, startingBehavior);
248     }
249     else if (observablesHistory->edsamHistory)
250     {
251         gmx_fatal(FARGS,
252                   "The checkpoint is from a run with essential dynamics sampling, "
253                   "but the current run did not specify the -ei option. "
254                   "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
255     }
256
257     initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
258     Update     upd(ir, deform);
259     const bool doSimulatedAnnealing = initSimulatedAnnealing(ir, &upd);
260     const bool useReplicaExchange   = (replExParams.exchangeInterval > 0);
261
262     bool simulationsShareState = false;
263     int  nstSignalComm         = nstglobalcomm;
264     {
265         // TODO This implementation of ensemble orientation restraints is nasty because
266         // a user can't just do multi-sim with single-sim orientation restraints.
267         bool usingEnsembleRestraints =
268                 (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
269         bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
270
271         // Replica exchange, ensemble restraints and AWH need all
272         // simulations to remain synchronized, so they need
273         // checkpoints and stop conditions to act on the same step, so
274         // the propagation of such signals must take place between
275         // simulations, not just within simulations.
276         // TODO: Make algorithm initializers set these flags.
277         simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
278
279         if (simulationsShareState)
280         {
281             // Inter-simulation signal communication does not need to happen
282             // often, so we use a minimum of 200 steps to reduce overhead.
283             const int c_minimumInterSimulationSignallingInterval = 200;
284             nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
285                             * nstglobalcomm;
286         }
287     }
288
289     if (startingBehavior != StartingBehavior::RestartWithAppending)
290     {
291         pleaseCiteCouplingAlgorithms(fplog, *ir);
292     }
293     gmx_mdoutf* outf =
294             init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, ir,
295                         top_global, oenv, wcycle, startingBehavior, simulationsShareState, ms);
296     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
297                                    mdoutf_get_fp_dhdl(outf), false, startingBehavior, mdModulesNotifier);
298
299     gstat = global_stat_init(ir);
300
301     /* Check for polarizable models and flexible constraints */
302     shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
303                                  ir->nstcalcenergy, DOMAINDECOMP(cr));
304
305     {
306         double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
307         if ((io > 2000) && MASTER(cr))
308         {
309             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
310         }
311     }
312
313     // Local state only becomes valid now.
314     std::unique_ptr<t_state> stateInstance;
315     t_state*                 state;
316
317
318     auto mdatoms = mdAtoms->mdatoms();
319
320     std::unique_ptr<UpdateConstrainCuda> integrator;
321
322     if (DOMAINDECOMP(cr))
323     {
324         dd_init_local_top(*top_global, &top);
325
326         stateInstance = std::make_unique<t_state>();
327         state         = stateInstance.get();
328         dd_init_local_state(cr->dd, state_global, state);
329
330         /* Distribute the charge groups over the nodes from the master node */
331         dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
332                             imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
333                             nrnb, nullptr, FALSE);
334         shouldCheckNumberOfBondedInteractions = true;
335         upd.setNumAtoms(state->natoms);
336     }
337     else
338     {
339         state_change_natoms(state_global, state_global->natoms);
340         f.resizeWithPadding(state_global->natoms);
341         /* Copy the pointer to the global state */
342         state = state_global;
343
344         /* Generate and initialize new topology */
345         mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &graph, mdAtoms, constr, vsite, shellfc);
346
347         upd.setNumAtoms(state->natoms);
348     }
349
350     const auto& simulationWork     = runScheduleWork->simulationWork;
351     const bool  useGpuForPme       = simulationWork.useGpuPme;
352     const bool  useGpuForNonbonded = simulationWork.useGpuNonbonded;
353     const bool  useGpuForBufferOps = simulationWork.useGpuBufferOps;
354     const bool  useGpuForUpdate    = simulationWork.useGpuUpdate;
355
356     StatePropagatorDataGpu* stateGpu = fr->stateGpu;
357
358     if (useGpuForUpdate)
359     {
360         GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
361                                    || constr->numConstraintsTotal() == 0,
362                            "Constraints in domain decomposition are only supported with update "
363                            "groups if using GPU update.\n");
364         GMX_RELEASE_ASSERT(ir->eConstrAlg != econtSHAKE || constr == nullptr
365                                    || constr->numConstraintsTotal() == 0,
366                            "SHAKE is not supported with GPU update.");
367         GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
368                            "Either PME or short-ranged non-bonded interaction tasks must run on "
369                            "the GPU to use GPU update.\n");
370         GMX_RELEASE_ASSERT(ir->eI == eiMD,
371                            "Only the md integrator is supported with the GPU update.\n");
372         GMX_RELEASE_ASSERT(
373                 ir->etc != etcNOSEHOOVER,
374                 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
375         GMX_RELEASE_ASSERT(ir->epc == epcNO || ir->epc == epcPARRINELLORAHMAN || ir->epc == epcBERENDSEN,
376                            "Only Parrinello-Rahman and Berendsen pressure coupling are supported "
377                            "with the GPU update.\n");
378         GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
379                            "Virtual sites are not supported with the GPU update.\n");
380         GMX_RELEASE_ASSERT(ed == nullptr,
381                            "Essential dynamics is not supported with the GPU update.\n");
382         GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(ir->pull),
383                            "Constraints pulling is not supported with the GPU update.\n");
384         GMX_RELEASE_ASSERT(fcd->orires.nr == 0,
385                            "Orientation restraints are not supported with the GPU update.\n");
386         GMX_RELEASE_ASSERT(ir->efep == efepNO,
387                            "Free energy perturbations are not supported with the GPU update.");
388         GMX_RELEASE_ASSERT(graph == nullptr, "The graph is not supported with GPU update.");
389
390         if (constr != nullptr && constr->numConstraintsTotal() > 0)
391         {
392             GMX_LOG(mdlog.info)
393                     .asParagraph()
394                     .appendText("Updating coordinates and applying constraints on the GPU.");
395         }
396         else
397         {
398             GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
399         }
400         integrator = std::make_unique<UpdateConstrainCuda>(
401                 *ir, *top_global, stateGpu->getUpdateStream(), stateGpu->xUpdatedOnDevice());
402
403         t_pbc pbc;
404         set_pbc(&pbc, epbcXYZ, state->box);
405         integrator->setPbc(&pbc);
406     }
407
408     if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
409     {
410         changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
411     }
412     if ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
413     {
414         changePinningPolicy(&f, PinningPolicy::PinnedIfSupported);
415     }
416     if (useGpuForUpdate)
417     {
418         changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
419     }
420
421     // NOTE: The global state is no longer used at this point.
422     // But state_global is still used as temporary storage space for writing
423     // the global state to file and potentially for replica exchange.
424     // (Global topology should persist.)
425
426     update_mdatoms(mdatoms, state->lambda[efptMASS]);
427
428     if (ir->bExpanded)
429     {
430         /* Check nstexpanded here, because the grompp check was broken */
431         if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
432         {
433             gmx_fatal(FARGS,
434                       "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
435         }
436         init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
437     }
438
439     if (MASTER(cr))
440     {
441         EnergyElement::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
442     }
443
444     preparePrevStepPullCom(ir, pull_work, mdatoms, state, state_global, cr,
445                            startingBehavior != StartingBehavior::NewSimulation);
446
447     // TODO: Remove this by converting AWH into a ForceProvider
448     auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms,
449                                 startingBehavior != StartingBehavior::NewSimulation,
450                                 shellfc != nullptr, opt2fn("-awh", nfile, fnm), pull_work);
451
452     if (useReplicaExchange && MASTER(cr))
453     {
454         repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir, replExParams);
455     }
456     /* PME tuning is only supported in the Verlet scheme, with PME for
457      * Coulomb. It is not supported with only LJ PME. */
458     bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
459                 && ir->cutoff_scheme != ecutsGROUP);
460
461     pme_load_balancing_t* pme_loadbal = nullptr;
462     if (bPMETune)
463     {
464         pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata,
465                          fr->nbv->useGpu());
466     }
467
468     if (!ir->bContinuation)
469     {
470         if (state->flags & (1U << estV))
471         {
472             auto v = makeArrayRef(state->v);
473             /* Set the velocities of vsites, shells and frozen atoms to zero */
474             for (i = 0; i < mdatoms->homenr; i++)
475             {
476                 if (mdatoms->ptype[i] == eptVSite || mdatoms->ptype[i] == eptShell)
477                 {
478                     clear_rvec(v[i]);
479                 }
480                 else if (mdatoms->cFREEZE)
481                 {
482                     for (m = 0; m < DIM; m++)
483                     {
484                         if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
485                         {
486                             v[i][m] = 0;
487                         }
488                     }
489                 }
490             }
491         }
492
493         if (constr)
494         {
495             /* Constrain the initial coordinates and velocities */
496             do_constrain_first(fplog, constr, ir, mdatoms, state->natoms, state->x.arrayRefWithPadding(),
497                                state->v.arrayRefWithPadding(), state->box, state->lambda[efptBONDED]);
498         }
499         if (vsite)
500         {
501             /* Construct the virtual sites for the initial configuration */
502             construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, nullptr, top.idef.iparams,
503                              top.idef.il, fr->ePBC, fr->bMolPBC, cr, state->box);
504         }
505     }
506
507     if (ir->efep != efepNO)
508     {
509         /* Set free energy calculation frequency as the greatest common
510          * denominator of nstdhdl and repl_ex_nst. */
511         nstfep = ir->fepvals->nstdhdl;
512         if (ir->bExpanded)
513         {
514             nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
515         }
516         if (useReplicaExchange)
517         {
518             nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
519         }
520     }
521
522     /* Be REALLY careful about what flags you set here. You CANNOT assume
523      * this is the first step, since we might be restarting from a checkpoint,
524      * and in that case we should not do any modifications to the state.
525      */
526     bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
527
528     // When restarting from a checkpoint, it can be appropriate to
529     // initialize ekind from quantities in the checkpoint. Otherwise,
530     // compute_globals must initialize ekind before the simulation
531     // starts/restarts. However, only the master rank knows what was
532     // found in the checkpoint file, so we have to communicate in
533     // order to coordinate the restart.
534     //
535     // TODO Consider removing this communication if/when checkpoint
536     // reading directly follows .tpr reading, because all ranks can
537     // agree on hasReadEkinState at that time.
538     bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
539     if (PAR(cr))
540     {
541         gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr);
542     }
543     if (hasReadEkinState)
544     {
545         restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
546     }
547
548     unsigned int cglo_flags =
549             (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
550              | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
551
552     bSumEkinhOld = FALSE;
553
554     t_vcm vcm(top_global->groups, *ir);
555     reportComRemovalInfo(fplog, vcm);
556
557     /* To minimize communication, compute_globals computes the COM velocity
558      * and the kinetic energy for the velocities without COM motion removed.
559      * Thus to get the kinetic energy without the COM contribution, we need
560      * to call compute_globals twice.
561      */
562     for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
563     {
564         unsigned int cglo_flags_iteration = cglo_flags;
565         if (bStopCM && cgloIteration == 0)
566         {
567             cglo_flags_iteration |= CGLO_STOPCM;
568             cglo_flags_iteration &= ~CGLO_TEMPERATURE;
569         }
570         compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
571                         state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, nullptr, enerd,
572                         force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
573                         state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
574                         cglo_flags_iteration
575                                 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
576                                                                          : 0));
577         if (cglo_flags_iteration & CGLO_STOPCM)
578         {
579             /* At initialization, do not pass x with acceleration-correction mode
580              * to avoid (incorrect) correction of the initial coordinates.
581              */
582             rvec* xPtr = nullptr;
583             if (vcm.mode != ecmLINEAR_ACCELERATION_CORRECTION)
584             {
585                 xPtr = state->x.rvec_array();
586             }
587             process_and_stopcm_grp(fplog, &vcm, *mdatoms, xPtr, state->v.rvec_array());
588             inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
589         }
590     }
591     checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
592                                     state->x.rvec_array(), state->box,
593                                     &shouldCheckNumberOfBondedInteractions);
594     if (ir->eI == eiVVAK)
595     {
596         /* a second call to get the half step temperature initialized as well */
597         /* we do the same call as above, but turn the pressure off -- internally to
598            compute_globals, this is recognized as a velocity verlet half-step
599            kinetic energy calculation.  This minimized excess variables, but
600            perhaps loses some logic?*/
601
602         compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
603                         state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, nullptr, enerd,
604                         force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
605                         state->box, nullptr, &bSumEkinhOld, cglo_flags & ~CGLO_PRESSURE);
606     }
607
608     /* Calculate the initial half step temperature, and save the ekinh_old */
609     if (startingBehavior == StartingBehavior::NewSimulation)
610     {
611         for (i = 0; (i < ir->opts.ngtc); i++)
612         {
613             copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
614         }
615     }
616
617     /* need to make an initiation call to get the Trotter variables set, as well as other constants
618        for non-trotter temperature control */
619     auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
620
621     if (MASTER(cr))
622     {
623         if (!ir->bContinuation)
624         {
625             if (constr && ir->eConstrAlg == econtLINCS)
626             {
627                 fprintf(fplog, "RMS relative constraint deviation after constraining: %.2e\n",
628                         constr->rmsd());
629             }
630             if (EI_STATE_VELOCITY(ir->eI))
631             {
632                 real temp = enerd->term[F_TEMP];
633                 if (ir->eI != eiVV)
634                 {
635                     /* Result of Ekin averaged over velocities of -half
636                      * and +half step, while we only have -half step here.
637                      */
638                     temp *= 2;
639                 }
640                 fprintf(fplog, "Initial temperature: %g K\n", temp);
641             }
642         }
643
644         char tbuf[20];
645         fprintf(stderr, "starting mdrun '%s'\n", *(top_global->name));
646         if (ir->nsteps >= 0)
647         {
648             sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
649         }
650         else
651         {
652             sprintf(tbuf, "%s", "infinite");
653         }
654         if (ir->init_step > 0)
655         {
656             fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
657                     gmx_step_str(ir->init_step + ir->nsteps, sbuf), tbuf,
658                     gmx_step_str(ir->init_step, sbuf2), ir->init_step * ir->delta_t);
659         }
660         else
661         {
662             fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
663         }
664         fprintf(fplog, "\n");
665     }
666
667     walltime_accounting_start_time(walltime_accounting);
668     wallcycle_start(wcycle, ewcRUN);
669     print_start(fplog, cr, walltime_accounting, "mdrun");
670
671 #if GMX_FAHCORE
672     /* safest point to do file checkpointing is here.  More general point would be immediately before integrator call */
673     int chkpt_ret = fcCheckPointParallel(cr->nodeid, NULL, 0);
674     if (chkpt_ret == 0)
675     {
676         gmx_fatal(3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0);
677     }
678 #endif
679
680     /***********************************************************
681      *
682      *             Loop over MD steps
683      *
684      ************************************************************/
685
686     bFirstStep = TRUE;
687     /* Skip the first Nose-Hoover integration when we get the state from tpx */
688     bInitStep        = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
689     bSumEkinhOld     = FALSE;
690     bExchanged       = FALSE;
691     bNeedRepartition = FALSE;
692
693     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
694             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
695             MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
696             mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
697
698     auto checkpointHandler = std::make_unique<CheckpointHandler>(
699             compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]), simulationsShareState,
700             ir->nstlist == 0, MASTER(cr), mdrunOptions.writeConfout,
701             mdrunOptions.checkpointOptions.period);
702
703     const bool resetCountersIsLocal = true;
704     auto       resetHandler         = std::make_unique<ResetHandler>(
705             compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
706             !resetCountersIsLocal, ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
707             mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
708
709     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
710
711     step     = ir->init_step;
712     step_rel = 0;
713
714     // TODO extract this to new multi-simulation module
715     if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
716     {
717         if (!multisim_int_all_are_equal(ms, ir->nsteps))
718         {
719             GMX_LOG(mdlog.warning)
720                     .appendText(
721                             "Note: The number of steps is not consistent across multi "
722                             "simulations,\n"
723                             "but we are proceeding anyway!");
724         }
725         if (!multisim_int_all_are_equal(ms, ir->init_step))
726         {
727             if (simulationsShareState)
728             {
729                 if (MASTER(cr))
730                 {
731                     gmx_fatal(FARGS,
732                               "The initial step is not consistent across multi simulations which "
733                               "share the state");
734                 }
735                 gmx_barrier(cr);
736             }
737             else
738             {
739                 GMX_LOG(mdlog.warning)
740                         .appendText(
741                                 "Note: The initial step is not consistent across multi "
742                                 "simulations,\n"
743                                 "but we are proceeding anyway!");
744             }
745         }
746     }
747
748     /* and stop now if we should */
749     bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
750     while (!bLastStep)
751     {
752
753         /* Determine if this is a neighbor search step */
754         bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
755
756         if (bPMETune && bNStList)
757         {
758             // This has to be here because PME load balancing is called so early.
759             // TODO: Move to after all booleans are defined.
760             if (useGpuForUpdate && !bFirstStep)
761             {
762                 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
763                 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
764             }
765             /* PME grid + cut-off optimization with GPUs or PME nodes */
766             pme_loadbal_do(pme_loadbal, cr, (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
767                            fplog, mdlog, *ir, fr, state->box, state->x, wcycle, step, step_rel,
768                            &bPMETunePrinting, simulationWork.useGpuPmePpCommunication);
769         }
770
771         wallcycle_start(wcycle, ewcSTEP);
772
773         bLastStep = (step_rel == ir->nsteps);
774         t         = t0 + step * ir->delta_t;
775
776         // TODO Refactor this, so that nstfep does not need a default value of zero
777         if (ir->efep != efepNO || ir->bSimTemp)
778         {
779             /* find and set the current lambdas */
780             setCurrentLambdasLocal(step, ir->fepvals, lam0, state->lambda, state->fep_state);
781
782             bDoDHDL     = do_per_step(step, ir->fepvals->nstdhdl);
783             bDoFEP      = ((ir->efep != efepNO) && do_per_step(step, nstfep));
784             bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
785                            && (step > 0) && (startingBehavior == StartingBehavior::NewSimulation));
786         }
787
788         bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
789                      && do_per_step(step, replExParams.exchangeInterval));
790
791         if (doSimulatedAnnealing)
792         {
793             update_annealing_target_temp(ir, t, &upd);
794         }
795
796         /* Stop Center of Mass motion */
797         bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
798
799         /* Determine whether or not to do Neighbour Searching */
800         bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
801
802         /* Note that the stopHandler will cause termination at nstglobalcomm
803          * steps. Since this concides with nstcalcenergy, nsttcouple and/or
804          * nstpcouple steps, we have computed the half-step kinetic energy
805          * of the previous step and can always output energies at the last step.
806          */
807         bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
808
809         /* do_log triggers energy and virial calculation. Because this leads
810          * to different code paths, forces can be different. Thus for exact
811          * continuation we should avoid extra log output.
812          * Note that the || bLastStep can result in non-exact continuation
813          * beyond the last step. But we don't consider that to be an issue.
814          */
815         do_log     = (do_per_step(step, ir->nstlog)
816                   || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
817         do_verbose = mdrunOptions.verbose
818                      && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
819
820         if (useGpuForUpdate && !bFirstStep && bNS)
821         {
822             // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
823             stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
824             stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
825             // Copy coordinate from the GPU when needed at the search step.
826             // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
827             // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
828             stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
829             stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
830         }
831
832         if (bNS && !(bFirstStep && ir->bContinuation))
833         {
834             bMasterState = FALSE;
835             /* Correct the new box if it is too skewed */
836             if (inputrecDynamicBox(ir))
837             {
838                 if (correct_box(fplog, step, state->box, graph))
839                 {
840                     bMasterState = TRUE;
841                     // If update is offloaded, it should be informed about the box size change
842                     if (useGpuForUpdate)
843                     {
844                         t_pbc pbc;
845                         set_pbc(&pbc, epbcXYZ, state->box);
846                         integrator->setPbc(&pbc);
847                     }
848                 }
849             }
850             if (DOMAINDECOMP(cr) && bMasterState)
851             {
852                 dd_collect_state(cr->dd, state, state_global);
853             }
854
855             if (DOMAINDECOMP(cr))
856             {
857                 /* Repartition the domain decomposition */
858                 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
859                                     *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
860                                     fr, vsite, constr, nrnb, wcycle, do_verbose && !bPMETunePrinting);
861                 shouldCheckNumberOfBondedInteractions = true;
862                 upd.setNumAtoms(state->natoms);
863             }
864         }
865
866         if (MASTER(cr) && do_log)
867         {
868             energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
869         }
870
871         if (ir->efep != efepNO)
872         {
873             update_mdatoms(mdatoms, state->lambda[efptMASS]);
874         }
875
876         if (bExchanged)
877         {
878
879             /* We need the kinetic energy at minus the half step for determining
880              * the full step kinetic energy and possibly for T-coupling.*/
881             /* This may not be quite working correctly yet . . . . */
882             compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
883                             state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
884                             nullptr, nullptr, nullptr, nullptr, mu_tot, constr, &nullSignaller,
885                             state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
886                             CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
887             checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
888                                             &top, state->x.rvec_array(), state->box,
889                                             &shouldCheckNumberOfBondedInteractions);
890         }
891         clear_mat(force_vir);
892
893         checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
894
895         /* Determine the energy and pressure:
896          * at nstcalcenergy steps and at energy output steps (set below).
897          */
898         if (EI_VV(ir->eI) && (!bInitStep))
899         {
900             bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
901             bCalcVir      = bCalcEnerStep
902                        || (ir->epc != epcNO
903                            && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
904         }
905         else
906         {
907             bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
908             bCalcVir = bCalcEnerStep || (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
909         }
910         bCalcEner = bCalcEnerStep;
911
912         do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
913
914         if (do_ene || do_log || bDoReplEx)
915         {
916             bCalcVir  = TRUE;
917             bCalcEner = TRUE;
918         }
919
920         /* Do we need global communication ? */
921         bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
922                   || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
923
924         force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
925                        | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
926                        | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
927
928         if (shellfc)
929         {
930             /* Now is the time to relax the shells */
931             relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
932                                 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
933                                 state->natoms, state->x.arrayRefWithPadding(),
934                                 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
935                                 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, graph,
936                                 shellfc, fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
937         }
938         else
939         {
940             /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
941                is updated (or the AWH update will be performed twice for one step when continuing).
942                It would be best to call this update function from do_md_trajectory_writing but that
943                would occur after do_force. One would have to divide the update_awh function into one
944                function applying the AWH force and one doing the AWH bias update. The update AWH
945                bias function could then be called after do_md_trajectory_writing (then containing
946                update_awh_history). The checkpointing will in the future probably moved to the start
947                of the md loop which will rid of this issue. */
948             if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
949             {
950                 awh->updateHistory(state_global->awhHistory.get());
951             }
952
953             /* The coordinates (x) are shifted (to get whole molecules)
954              * in do_force.
955              * This is parallellized as well, and does communication too.
956              * Check comments in sim_util.c
957              */
958             do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation, imdSession, pull_work, step,
959                      nrnb, wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
960                      f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, graph,
961                      fr, runScheduleWork, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
962                      (bNS ? GMX_FORCE_NS : 0) | force_flags, ddBalanceRegionHandler);
963         }
964
965         // VV integrators do not need the following velocity half step
966         // if it is the first step after starting from a checkpoint.
967         // That is, the half step is needed on all other steps, and
968         // also the first step when starting from a .tpr file.
969         if (EI_VV(ir->eI) && (!bFirstStep || startingBehavior == StartingBehavior::NewSimulation))
970         /*  ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
971         {
972             rvec* vbuf = nullptr;
973
974             wallcycle_start(wcycle, ewcUPDATE);
975             if (ir->eI == eiVV && bInitStep)
976             {
977                 /* if using velocity verlet with full time step Ekin,
978                  * take the first half step only to compute the
979                  * virial for the first step. From there,
980                  * revert back to the initial coordinates
981                  * so that the input is actually the initial step.
982                  */
983                 snew(vbuf, state->natoms);
984                 copy_rvecn(state->v.rvec_array(), vbuf, 0,
985                            state->natoms); /* should make this better for parallelizing? */
986             }
987             else
988             {
989                 /* this is for NHC in the Ekin(t+dt/2) version of vv */
990                 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
991                                trotter_seq, ettTSEQ1);
992             }
993
994             update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
995                           etrtVELOCITY1, cr, constr);
996
997             wallcycle_stop(wcycle, ewcUPDATE);
998             constrain_velocities(step, nullptr, state, shake_vir, constr, bCalcVir, do_log, do_ene);
999             wallcycle_start(wcycle, ewcUPDATE);
1000             /* if VV, compute the pressure and constraints */
1001             /* For VV2, we strictly only need this if using pressure
1002              * control, but we really would like to have accurate pressures
1003              * printed out.
1004              * Think about ways around this in the future?
1005              * For now, keep this choice in comments.
1006              */
1007             /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1008             /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1009             bPres = TRUE;
1010             bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1011             if (bCalcEner && ir->eI == eiVVAK)
1012             {
1013                 bSumEkinhOld = TRUE;
1014             }
1015             /* for vv, the first half of the integration actually corresponds to the previous step.
1016                So we need information from the last step in the first half of the integration */
1017             if (bGStat || do_per_step(step - 1, nstglobalcomm))
1018             {
1019                 wallcycle_stop(wcycle, ewcUPDATE);
1020                 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1021                                 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1022                                 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
1023                                 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1024                                 (bGStat ? CGLO_GSTAT : 0) | (bCalcEner ? CGLO_ENERGY : 0)
1025                                         | (bTemp ? CGLO_TEMPERATURE : 0) | (bPres ? CGLO_PRESSURE : 0)
1026                                         | (bPres ? CGLO_CONSTRAINT : 0) | (bStopCM ? CGLO_STOPCM : 0)
1027                                         | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1028                                                                                  : 0)
1029                                         | CGLO_SCALEEKIN);
1030                 /* explanation of above:
1031                    a) We compute Ekin at the full time step
1032                    if 1) we are using the AveVel Ekin, and it's not the
1033                    initial step, or 2) if we are using AveEkin, but need the full
1034                    time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1035                    b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1036                    EkinAveVel because it's needed for the pressure */
1037                 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1038                                                 top_global, &top, state->x.rvec_array(), state->box,
1039                                                 &shouldCheckNumberOfBondedInteractions);
1040                 if (bStopCM)
1041                 {
1042                     process_and_stopcm_grp(fplog, &vcm, *mdatoms, state->x.rvec_array(),
1043                                            state->v.rvec_array());
1044                     inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1045                 }
1046                 wallcycle_start(wcycle, ewcUPDATE);
1047             }
1048             /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1049             if (!bInitStep)
1050             {
1051                 if (bTrotter)
1052                 {
1053                     m_add(force_vir, shake_vir,
1054                           total_vir); /* we need the un-dispersion corrected total vir here */
1055                     trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1056                                    trotter_seq, ettTSEQ2);
1057
1058                     /* TODO This is only needed when we're about to write
1059                      * a checkpoint, because we use it after the restart
1060                      * (in a kludge?). But what should we be doing if
1061                      * the startingBehavior is NewSimulation or bInitStep are true? */
1062                     if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1063                     {
1064                         copy_mat(shake_vir, state->svir_prev);
1065                         copy_mat(force_vir, state->fvir_prev);
1066                     }
1067                     if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1068                     {
1069                         /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1070                         enerd->term[F_TEMP] =
1071                                 sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1072                         enerd->term[F_EKIN] = trace(ekind->ekin);
1073                     }
1074                 }
1075                 else if (bExchanged)
1076                 {
1077                     wallcycle_stop(wcycle, ewcUPDATE);
1078                     /* We need the kinetic energy at minus the half step for determining
1079                      * the full step kinetic energy and possibly for T-coupling.*/
1080                     /* This may not be quite working correctly yet . . . . */
1081                     compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(),
1082                                     state->v.rvec_array(), state->box, state->lambda[efptVDW],
1083                                     mdatoms, nrnb, &vcm, wcycle, enerd, nullptr, nullptr, nullptr,
1084                                     nullptr, mu_tot, constr, &nullSignaller, state->box, nullptr,
1085                                     &bSumEkinhOld, CGLO_GSTAT | CGLO_TEMPERATURE);
1086                     wallcycle_start(wcycle, ewcUPDATE);
1087                 }
1088             }
1089             /* if it's the initial step, we performed this first step just to get the constraint virial */
1090             if (ir->eI == eiVV && bInitStep)
1091             {
1092                 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1093                 sfree(vbuf);
1094             }
1095             wallcycle_stop(wcycle, ewcUPDATE);
1096         }
1097
1098         /* compute the conserved quantity */
1099         if (EI_VV(ir->eI))
1100         {
1101             saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1102             if (ir->eI == eiVV)
1103             {
1104                 last_ekin = enerd->term[F_EKIN];
1105             }
1106             if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1107             {
1108                 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1109             }
1110             /* sum up the foreign energy and dhdl terms for vv.  currently done every step so that dhdl is correct in the .edr */
1111             if (ir->efep != efepNO)
1112             {
1113                 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1114             }
1115         }
1116
1117         /* ########  END FIRST UPDATE STEP  ############## */
1118         /* ########  If doing VV, we now have v(dt) ###### */
1119         if (bDoExpanded)
1120         {
1121             /* perform extended ensemble sampling in lambda - we don't
1122                actually move to the new state before outputting
1123                statistics, but if performing simulated tempering, we
1124                do update the velocities and the tau_t. */
1125
1126             lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state,
1127                                               state->dfhist, step, state->v.rvec_array(), mdatoms);
1128             /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1129             if (MASTER(cr))
1130             {
1131                 copy_df_history(state_global->dfhist, state->dfhist);
1132             }
1133         }
1134
1135         // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1136         // coordinates have not already been copied for i) search or ii) CPU force tasks.
1137         if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1138             && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1139                 || checkpointHandler->isCheckpointingStep()))
1140         {
1141             stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1142             stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1143         }
1144         // Copy velocities if needed for the output/checkpointing.
1145         // NOTE: Copy on the search steps is done at the beginning of the step.
1146         if (useGpuForUpdate && !bNS
1147             && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1148         {
1149             stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1150             stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1151         }
1152         // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1153         // and update is offloaded hence forces are kept on the GPU for update and have not been
1154         // already transferred in do_force().
1155         // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1156         //       when the forces are ready on the GPU -- the same synchronizer should be used as the one
1157         //       prior to GPU update.
1158         // TODO: When the output flags will be included in step workload, this copy can be combined with the
1159         //       copy call in do_force(...).
1160         // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1161         //       on host after the D2H copy in do_force(...).
1162         if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1163             && do_per_step(step, ir->nstfout))
1164         {
1165             stateGpu->copyForcesFromGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1166             stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1167         }
1168         /* Now we have the energies and forces corresponding to the
1169          * coordinates at time t. We must output all of this before
1170          * the update.
1171          */
1172         do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state, state_global,
1173                                  observablesHistory, top_global, fr, outf, energyOutput, ekind, f,
1174                                  checkpointHandler->isCheckpointingStep(), bRerunMD, bLastStep,
1175                                  mdrunOptions.writeConfout, bSumEkinhOld);
1176         /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1177         bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1178
1179         /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1180         if (startingBehavior != StartingBehavior::NewSimulation
1181             && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1182         {
1183             copy_mat(state->svir_prev, shake_vir);
1184             copy_mat(state->fvir_prev, force_vir);
1185         }
1186
1187         stopHandler->setSignal();
1188         resetHandler->setSignal(walltime_accounting);
1189
1190         if (bGStat || !PAR(cr))
1191         {
1192             /* In parallel we only have to check for checkpointing in steps
1193              * where we do global communication,
1194              *  otherwise the other nodes don't know.
1195              */
1196             checkpointHandler->setSignal(walltime_accounting);
1197         }
1198
1199         /* #########   START SECOND UPDATE STEP ################# */
1200
1201         /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1202            controlled in preprocessing */
1203
1204         if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1205         {
1206             gmx_bool bIfRandomize;
1207             bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1208             /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1209             if (constr && bIfRandomize)
1210             {
1211                 constrain_velocities(step, nullptr, state, tmp_vir, constr, bCalcVir, do_log, do_ene);
1212             }
1213         }
1214         /* Box is changed in update() when we do pressure coupling,
1215          * but we should still use the old box for energy corrections and when
1216          * writing it to the energy file, so it matches the trajectory files for
1217          * the same timestep above. Make a copy in a separate array.
1218          */
1219         copy_mat(state->box, lastbox);
1220
1221         dvdl_constr = 0;
1222
1223         wallcycle_start(wcycle, ewcUPDATE);
1224         /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1225         if (bTrotter)
1226         {
1227             trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1228             /* We can only do Berendsen coupling after we have summed
1229              * the kinetic energy or virial. Since the happens
1230              * in global_state after update, we should only do it at
1231              * step % nstlist = 1 with bGStatEveryStep=FALSE.
1232              */
1233         }
1234         else
1235         {
1236             update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1237             update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1238         }
1239
1240         if (EI_VV(ir->eI))
1241         {
1242             /* velocity half-step update */
1243             update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1244                           etrtVELOCITY2, cr, constr);
1245         }
1246
1247         /* Above, initialize just copies ekinh into ekin,
1248          * it doesn't copy position (for VV),
1249          * and entire integrator for MD.
1250          */
1251
1252         if (ir->eI == eiVVAK)
1253         {
1254             cbuf.resize(state->x.size());
1255             std::copy(state->x.begin(), state->x.end(), cbuf.begin());
1256         }
1257
1258         /* With leap-frog type integrators we compute the kinetic energy
1259          * at a whole time step as the average of the half-time step kinetic
1260          * energies of two subsequent steps. Therefore we need to compute the
1261          * half step kinetic energy also if we need energies at the next step.
1262          */
1263         const bool needHalfStepKineticEnergy =
1264                 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1265
1266         // Parrinello-Rahman requires the pressure to be availible before the update to compute
1267         // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1268         const bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN
1269                                          && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1270
1271         if (useGpuForUpdate)
1272         {
1273             if (bNS && (bFirstStep || DOMAINDECOMP(cr)))
1274             {
1275                 integrator->set(stateGpu->getCoordinates(), stateGpu->getVelocities(),
1276                                 stateGpu->getForces(), top.idef, *mdatoms, ekind->ngtc);
1277
1278                 // Copy data to the GPU after buffers might have being reinitialized
1279                 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1280                 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1281             }
1282
1283             // If the buffer ops were not offloaded this step, the forces are on the host and have to be copied
1284             if (!runScheduleWork->stepWork.useGpuFBufferOps)
1285             {
1286                 stateGpu->copyForcesToGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1287             }
1288
1289             const bool doTemperatureScaling =
1290                     (ir->etc != etcNO && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1291
1292             // This applies Leap-Frog, LINCS and SETTLE in succession
1293             integrator->integrate(stateGpu->getForcesReadyOnDeviceEvent(
1294                                           AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1295                                   ir->delta_t, true, bCalcVir, shake_vir, doTemperatureScaling,
1296                                   ekind->tcstat, doParrinelloRahman, ir->nstpcouple * ir->delta_t, M);
1297
1298             // Copy velocities D2H after update if:
1299             // - Globals are computed this step (includes the energy output steps).
1300             // - Temperature is needed for the next step.
1301             if (bGStat || needHalfStepKineticEnergy)
1302             {
1303                 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1304                 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1305             }
1306         }
1307         else
1308         {
1309             update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1310                           etrtPOSITION, cr, constr);
1311
1312             wallcycle_stop(wcycle, ewcUPDATE);
1313
1314             constrain_coordinates(step, &dvdl_constr, state, shake_vir, &upd, constr, bCalcVir,
1315                                   do_log, do_ene);
1316
1317             update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state, cr, nrnb, wcycle, &upd,
1318                                   constr, do_log, do_ene);
1319             finish_update(ir, mdatoms, state, graph, nrnb, wcycle, &upd, constr);
1320         }
1321
1322         if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1323         {
1324             updatePrevStepPullCom(pull_work, state);
1325         }
1326
1327         if (ir->eI == eiVVAK)
1328         {
1329             /* erase F_EKIN and F_TEMP here? */
1330             /* just compute the kinetic energy at the half step to perform a trotter step */
1331             compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1332                             state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1333                             force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller, lastbox,
1334                             nullptr, &bSumEkinhOld, (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE);
1335             wallcycle_start(wcycle, ewcUPDATE);
1336             trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1337             /* now we know the scaling, we can compute the positions again */
1338             std::copy(cbuf.begin(), cbuf.end(), state->x.begin());
1339
1340             update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1341                           etrtPOSITION, cr, constr);
1342             wallcycle_stop(wcycle, ewcUPDATE);
1343
1344             /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1345             /* are the small terms in the shake_vir here due
1346              * to numerical errors, or are they important
1347              * physically? I'm thinking they are just errors, but not completely sure.
1348              * For now, will call without actually constraining, constr=NULL*/
1349             finish_update(ir, mdatoms, state, graph, nrnb, wcycle, &upd, nullptr);
1350         }
1351         if (EI_VV(ir->eI))
1352         {
1353             /* this factor or 2 correction is necessary
1354                because half of the constraint force is removed
1355                in the vv step, so we have to double it.  See
1356                the Redmine issue #1255.  It is not yet clear
1357                if the factor of 2 is exact, or just a very
1358                good approximation, and this will be
1359                investigated.  The next step is to see if this
1360                can be done adding a dhdl contribution from the
1361                rattle step, but this is somewhat more
1362                complicated with the current code. Will be
1363                investigated, hopefully for 4.6.3. However,
1364                this current solution is much better than
1365                having it completely wrong.
1366              */
1367             enerd->term[F_DVDL_CONSTR] += 2 * dvdl_constr;
1368         }
1369         else
1370         {
1371             enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1372         }
1373
1374         if (vsite != nullptr)
1375         {
1376             wallcycle_start(wcycle, ewcVSITECONSTR);
1377             if (graph != nullptr)
1378             {
1379                 shift_self(graph, state->box, state->x.rvec_array());
1380             }
1381             construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, state->v.rvec_array(),
1382                              top.idef.iparams, top.idef.il, fr->ePBC, fr->bMolPBC, cr, state->box);
1383
1384             if (graph != nullptr)
1385             {
1386                 unshift_self(graph, state->box, state->x.rvec_array());
1387             }
1388             wallcycle_stop(wcycle, ewcVSITECONSTR);
1389         }
1390
1391         /* ############## IF NOT VV, Calculate globals HERE  ############ */
1392         /* With Leap-Frog we can skip compute_globals at
1393          * non-communication steps, but we need to calculate
1394          * the kinetic energy one step before communication.
1395          */
1396         {
1397             // Organize to do inter-simulation signalling on steps if
1398             // and when algorithms require it.
1399             const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1400
1401             if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1402             {
1403                 // Copy coordinates when needed to stop the CM motion.
1404                 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1405                 {
1406                     stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1407                     stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1408                 }
1409                 // Since we're already communicating at this step, we
1410                 // can propagate intra-simulation signals. Note that
1411                 // check_nstglobalcomm has the responsibility for
1412                 // choosing the value of nstglobalcomm that is one way
1413                 // bGStat becomes true, so we can't get into a
1414                 // situation where e.g. checkpointing can't be
1415                 // signalled.
1416                 bool                doIntraSimSignal = true;
1417                 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1418
1419                 compute_globals(
1420                         gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1421                         state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1422                         force_vir, shake_vir, total_vir, pres, mu_tot, constr, &signaller, lastbox,
1423                         &totalNumberOfBondedInteractions, &bSumEkinhOld,
1424                         (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1425                                 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1426                                 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1427                                 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1428                                 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1429                                                                          : 0));
1430                 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1431                                                 top_global, &top, state->x.rvec_array(), state->box,
1432                                                 &shouldCheckNumberOfBondedInteractions);
1433                 if (!EI_VV(ir->eI) && bStopCM)
1434                 {
1435                     process_and_stopcm_grp(fplog, &vcm, *mdatoms, state->x.rvec_array(),
1436                                            state->v.rvec_array());
1437                     inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1438
1439                     // TODO: The special case of removing CM motion should be dealt more gracefully
1440                     if (useGpuForUpdate)
1441                     {
1442                         stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1443                         // Here we block until the H2D copy completes because event sync with the
1444                         // force kernels that use the coordinates on the next steps is not implemented
1445                         // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1446                         stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1447                     }
1448                 }
1449             }
1450         }
1451
1452         /* #############  END CALC EKIN AND PRESSURE ################# */
1453
1454         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1455            the virial that should probably be addressed eventually. state->veta has better properies,
1456            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1457            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
1458
1459         if (ir->efep != efepNO && !EI_VV(ir->eI))
1460         {
1461             /* Sum up the foreign energy and dhdl terms for md and sd.
1462                Currently done every step so that dhdl is correct in the .edr */
1463             sum_dhdl(enerd, state->lambda, *ir->fepvals);
1464         }
1465
1466         update_pcouple_after_coordinates(fplog, step, ir, mdatoms, pres, force_vir, shake_vir,
1467                                          pressureCouplingMu, state, nrnb, &upd, !useGpuForUpdate);
1468
1469         const bool doBerendsenPressureCoupling =
1470                 (inputrec->epc == epcBERENDSEN && do_per_step(step, inputrec->nstpcouple));
1471         if (useGpuForUpdate && (doBerendsenPressureCoupling || doParrinelloRahman))
1472         {
1473             integrator->scaleCoordinates(pressureCouplingMu);
1474             t_pbc pbc;
1475             set_pbc(&pbc, epbcXYZ, state->box);
1476             integrator->setPbc(&pbc);
1477         }
1478
1479         /* ################# END UPDATE STEP 2 ################# */
1480         /* #### We now have r(t+dt) and v(t+dt/2)  ############# */
1481
1482         /* The coordinates (x) were unshifted in update */
1483         if (!bGStat)
1484         {
1485             /* We will not sum ekinh_old,
1486              * so signal that we still have to do it.
1487              */
1488             bSumEkinhOld = TRUE;
1489         }
1490
1491         if (bCalcEner)
1492         {
1493             /* #########  BEGIN PREPARING EDR OUTPUT  ###########  */
1494
1495             /* use the directly determined last velocity, not actually the averaged half steps */
1496             if (bTrotter && ir->eI == eiVV)
1497             {
1498                 enerd->term[F_EKIN] = last_ekin;
1499             }
1500             enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1501
1502             if (integratorHasConservedEnergyQuantity(ir))
1503             {
1504                 if (EI_VV(ir->eI))
1505                 {
1506                     enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1507                 }
1508                 else
1509                 {
1510                     enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1511                 }
1512             }
1513             /* #########  END PREPARING EDR OUTPUT  ###########  */
1514         }
1515
1516         /* Output stuff */
1517         if (MASTER(cr))
1518         {
1519             if (fplog && do_log && bDoExpanded)
1520             {
1521                 /* only needed if doing expanded ensemble */
1522                 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals,
1523                                           ir->bSimTemp ? ir->simtempvals : nullptr,
1524                                           state_global->dfhist, state->fep_state, ir->nstlog, step);
1525             }
1526             if (bCalcEner)
1527             {
1528                 energyOutput.addDataAtEnergyStep(bDoDHDL, bCalcEnerStep, t, mdatoms->tmass, enerd, state,
1529                                                  ir->fepvals, ir->expandedvals, lastbox, shake_vir,
1530                                                  force_vir, total_vir, pres, ekind, mu_tot, constr);
1531             }
1532             else
1533             {
1534                 energyOutput.recordNonEnergyStep();
1535             }
1536
1537             gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1538             gmx_bool do_or = do_per_step(step, ir->nstorireout);
1539
1540             if (doSimulatedAnnealing)
1541             {
1542                 energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
1543             }
1544             if (do_log || do_ene || do_dr || do_or)
1545             {
1546                 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1547                                                    do_log ? fplog : nullptr, step, t, fcd, awh.get());
1548             }
1549
1550             if (ir->bPull)
1551             {
1552                 pull_print_output(pull_work, step, t);
1553             }
1554
1555             if (do_per_step(step, ir->nstlog))
1556             {
1557                 if (fflush(fplog) != 0)
1558                 {
1559                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1560                 }
1561             }
1562         }
1563         if (bDoExpanded)
1564         {
1565             /* Have to do this part _after_ outputting the logfile and the edr file */
1566             /* Gets written into the state at the beginning of next loop*/
1567             state->fep_state = lamnew;
1568         }
1569         /* Print the remaining wall clock time for the run */
1570         if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1571         {
1572             if (shellfc)
1573             {
1574                 fprintf(stderr, "\n");
1575             }
1576             print_time(stderr, walltime_accounting, step, ir, cr);
1577         }
1578
1579         /* Ion/water position swapping.
1580          * Not done in last step since trajectory writing happens before this call
1581          * in the MD loop and exchanges would be lost anyway. */
1582         bNeedRepartition = FALSE;
1583         if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep && do_per_step(step, ir->swap->nstswap))
1584         {
1585             bNeedRepartition =
1586                     do_swapcoords(cr, step, t, ir, swap, wcycle, as_rvec_array(state->x.data()),
1587                                   state->box, MASTER(cr) && mdrunOptions.verbose, bRerunMD);
1588
1589             if (bNeedRepartition && DOMAINDECOMP(cr))
1590             {
1591                 dd_collect_state(cr->dd, state, state_global);
1592             }
1593         }
1594
1595         /* Replica exchange */
1596         bExchanged = FALSE;
1597         if (bDoReplEx)
1598         {
1599             bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1600         }
1601
1602         if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1603         {
1604             dd_partition_system(fplog, mdlog, step, cr, TRUE, 1, state_global, *top_global, ir,
1605                                 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
1606                                 nrnb, wcycle, FALSE);
1607             shouldCheckNumberOfBondedInteractions = true;
1608             upd.setNumAtoms(state->natoms);
1609         }
1610
1611         bFirstStep = FALSE;
1612         bInitStep  = FALSE;
1613
1614         /* #######  SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1615         /* With all integrators, except VV, we need to retain the pressure
1616          * at the current step for coupling at the next step.
1617          */
1618         if ((state->flags & (1U << estPRES_PREV))
1619             && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1620         {
1621             /* Store the pressure in t_state for pressure coupling
1622              * at the next MD step.
1623              */
1624             copy_mat(pres, state->pres_prev);
1625         }
1626
1627         /* #######  END SET VARIABLES FOR NEXT ITERATION ###### */
1628
1629         if ((membed != nullptr) && (!bLastStep))
1630         {
1631             rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1632         }
1633
1634         cycles = wallcycle_stop(wcycle, ewcSTEP);
1635         if (DOMAINDECOMP(cr) && wcycle)
1636         {
1637             dd_cycles_add(cr->dd, cycles, ddCyclStep);
1638         }
1639
1640         /* increase the MD step number */
1641         step++;
1642         step_rel++;
1643
1644         resetHandler->resetCounters(step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb,
1645                                     fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1646
1647         /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1648         imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1649     }
1650     /* End of main MD loop */
1651
1652     /* Closing TNG files can include compressing data. Therefore it is good to do that
1653      * before stopping the time measurements. */
1654     mdoutf_tng_close(outf);
1655
1656     /* Stop measuring walltime */
1657     walltime_accounting_end_time(walltime_accounting);
1658
1659     if (!thisRankHasDuty(cr, DUTY_PME))
1660     {
1661         /* Tell the PME only node to finish */
1662         gmx_pme_send_finish(cr);
1663     }
1664
1665     if (MASTER(cr))
1666     {
1667         if (ir->nstcalcenergy > 0)
1668         {
1669             energyOutput.printAnnealingTemperatures(fplog, groups, &(ir->opts));
1670             energyOutput.printAverages(fplog, groups);
1671         }
1672     }
1673     done_mdoutf(outf);
1674
1675     if (bPMETune)
1676     {
1677         pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1678     }
1679
1680     done_shellfc(fplog, shellfc, step_rel);
1681
1682     if (useReplicaExchange && MASTER(cr))
1683     {
1684         print_replica_exchange_statistics(fplog, repl_ex);
1685     }
1686
1687     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1688
1689     global_stat_destroy(gstat);
1690 }