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