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