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