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