Merge branch release-5-1 into release-2016
[alexxy/gromacs.git] / src / programs / 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,2012,2013,2014,2015,2016, 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 #include "gmxpre.h"
38
39 #include "md.h"
40
41 #include "config.h"
42
43 #include <math.h>
44 #include <stdio.h>
45 #include <stdlib.h>
46
47 #include <algorithm>
48
49 #include "thread_mpi/threads.h"
50
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_network.h"
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/ewald/pme.h"
56 #include "gromacs/ewald/pme-load-balancing.h"
57 #include "gromacs/fileio/trxio.h"
58 #include "gromacs/gmxlib/md_logging.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/gmxlib/nrnb.h"
61 #include "gromacs/gpu_utils/gpu_utils.h"
62 #include "gromacs/imd/imd.h"
63 #include "gromacs/listed-forces/manage-threading.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/utilities.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/math/vectypes.h"
68 #include "gromacs/mdlib/compute_io.h"
69 #include "gromacs/mdlib/constr.h"
70 #include "gromacs/mdlib/ebin.h"
71 #include "gromacs/mdlib/force.h"
72 #include "gromacs/mdlib/force_flags.h"
73 #include "gromacs/mdlib/forcerec.h"
74 #include "gromacs/mdlib/md_support.h"
75 #include "gromacs/mdlib/mdatoms.h"
76 #include "gromacs/mdlib/mdebin.h"
77 #include "gromacs/mdlib/mdoutf.h"
78 #include "gromacs/mdlib/mdrun.h"
79 #include "gromacs/mdlib/nb_verlet.h"
80 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
81 #include "gromacs/mdlib/ns.h"
82 #include "gromacs/mdlib/shellfc.h"
83 #include "gromacs/mdlib/sighandler.h"
84 #include "gromacs/mdlib/sim_util.h"
85 #include "gromacs/mdlib/simulationsignal.h"
86 #include "gromacs/mdlib/tgroup.h"
87 #include "gromacs/mdlib/trajectory_writing.h"
88 #include "gromacs/mdlib/update.h"
89 #include "gromacs/mdlib/vcm.h"
90 #include "gromacs/mdlib/vsite.h"
91 #include "gromacs/mdtypes/commrec.h"
92 #include "gromacs/mdtypes/df_history.h"
93 #include "gromacs/mdtypes/energyhistory.h"
94 #include "gromacs/mdtypes/fcdata.h"
95 #include "gromacs/mdtypes/forcerec.h"
96 #include "gromacs/mdtypes/group.h"
97 #include "gromacs/mdtypes/inputrec.h"
98 #include "gromacs/mdtypes/interaction_const.h"
99 #include "gromacs/mdtypes/md_enums.h"
100 #include "gromacs/mdtypes/mdatom.h"
101 #include "gromacs/mdtypes/state.h"
102 #include "gromacs/pbcutil/mshift.h"
103 #include "gromacs/pbcutil/pbc.h"
104 #include "gromacs/pulling/pull.h"
105 #include "gromacs/swap/swapcoords.h"
106 #include "gromacs/timing/wallcycle.h"
107 #include "gromacs/timing/walltime_accounting.h"
108 #include "gromacs/topology/atoms.h"
109 #include "gromacs/topology/idef.h"
110 #include "gromacs/topology/mtop_util.h"
111 #include "gromacs/topology/topology.h"
112 #include "gromacs/trajectory/trajectoryframe.h"
113 #include "gromacs/utility/basedefinitions.h"
114 #include "gromacs/utility/cstringutil.h"
115 #include "gromacs/utility/fatalerror.h"
116 #include "gromacs/utility/real.h"
117 #include "gromacs/utility/smalloc.h"
118
119 #include "deform.h"
120 #include "membed.h"
121 #include "repl_ex.h"
122
123 #ifdef GMX_FAHCORE
124 #include "corewrap.h"
125 #endif
126
127 using gmx::SimulationSignaller;
128
129 /*! \brief Check whether bonded interactions are missing, if appropriate
130  *
131  * \param[in]    fplog                                  Log file pointer
132  * \param[in]    cr                                     Communication object
133  * \param[in]    totalNumberOfBondedInteractions        Result of the global reduction over the number of bonds treated in each domain
134  * \param[in]    top_global                             Global topology for the error message
135  * \param[in]    top_local                              Local topology for the error message
136  * \param[in]    state                                  Global state for the error message
137  * \param[inout] shouldCheckNumberOfBondedInteractions  Whether we should do the check.
138  *
139  * \return Nothing, except that shouldCheckNumberOfBondedInteractions
140  * is always set to false after exit.
141  */
142 static void checkNumberOfBondedInteractions(FILE *fplog, t_commrec *cr, int totalNumberOfBondedInteractions,
143                                             gmx_mtop_t *top_global, gmx_localtop_t *top_local, t_state *state,
144                                             bool *shouldCheckNumberOfBondedInteractions)
145 {
146     if (*shouldCheckNumberOfBondedInteractions)
147     {
148         if (totalNumberOfBondedInteractions != cr->dd->nbonded_global)
149         {
150             dd_print_missing_interactions(fplog, cr, totalNumberOfBondedInteractions, top_global, top_local, state); // Does not return
151         }
152         *shouldCheckNumberOfBondedInteractions = false;
153     }
154 }
155
156 static void reset_all_counters(FILE *fplog, t_commrec *cr,
157                                gmx_int64_t step,
158                                gmx_int64_t *step_rel, t_inputrec *ir,
159                                gmx_wallcycle_t wcycle, t_nrnb *nrnb,
160                                gmx_walltime_accounting_t walltime_accounting,
161                                struct nonbonded_verlet_t *nbv)
162 {
163     char sbuf[STEPSTRSIZE];
164
165     /* Reset all the counters related to performance over the run */
166     md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
167                   gmx_step_str(step, sbuf));
168
169     if (use_GPU(nbv))
170     {
171         nbnxn_gpu_reset_timings(nbv);
172         resetGpuProfiler();
173     }
174
175     wallcycle_stop(wcycle, ewcRUN);
176     wallcycle_reset_all(wcycle);
177     if (DOMAINDECOMP(cr))
178     {
179         reset_dd_statistics_counters(cr->dd);
180     }
181     init_nrnb(nrnb);
182     ir->init_step += *step_rel;
183     ir->nsteps    -= *step_rel;
184     *step_rel      = 0;
185     wallcycle_start(wcycle, ewcRUN);
186     walltime_accounting_start(walltime_accounting);
187     print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
188 }
189
190 /*! \libinternal
191     \copydoc integrator_t (FILE *fplog, t_commrec *cr,
192                            int nfile, const t_filenm fnm[],
193                            const gmx_output_env_t *oenv, gmx_bool bVerbose,
194                            int nstglobalcomm,
195                            gmx_vsite_t *vsite, gmx_constr_t constr,
196                            int stepout,
197                            t_inputrec *inputrec,
198                            gmx_mtop_t *top_global, t_fcdata *fcd,
199                            t_state *state_global,
200                            t_mdatoms *mdatoms,
201                            t_nrnb *nrnb, gmx_wallcycle_t wcycle,
202                            gmx_edsam_t ed,
203                            t_forcerec *fr,
204                            int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
205                            real cpt_period, real max_hours,
206                            int imdport,
207                            unsigned long Flags,
208                            gmx_walltime_accounting_t walltime_accounting)
209  */
210 double gmx::do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
211                   const gmx_output_env_t *oenv, gmx_bool bVerbose,
212                   int nstglobalcomm,
213                   gmx_vsite_t *vsite, gmx_constr_t constr,
214                   int stepout, t_inputrec *ir,
215                   gmx_mtop_t *top_global,
216                   t_fcdata *fcd,
217                   t_state *state_global,
218                   t_mdatoms *mdatoms,
219                   t_nrnb *nrnb, gmx_wallcycle_t wcycle,
220                   gmx_edsam_t ed, t_forcerec *fr,
221                   int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
222                   gmx_membed_t *membed,
223                   real cpt_period, real max_hours,
224                   int imdport,
225                   unsigned long Flags,
226                   gmx_walltime_accounting_t walltime_accounting)
227 {
228     gmx_mdoutf_t    outf = NULL;
229     gmx_int64_t     step, step_rel;
230     double          elapsed_time;
231     double          t, t0, lam0[efptNR];
232     gmx_bool        bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
233     gmx_bool        bNS, bNStList, bSimAnn, bStopCM, bRerunMD,
234                     bFirstStep, startingFromCheckpoint, bInitStep, bLastStep = FALSE,
235                     bBornRadii, bUsingEnsembleRestraints;
236     gmx_bool          bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
237     gmx_bool          do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
238                       bForceUpdate = FALSE, bCPT;
239     gmx_bool          bMasterState;
240     int               force_flags, cglo_flags;
241     tensor            force_vir, shake_vir, total_vir, tmp_vir, pres;
242     int               i, m;
243     t_trxstatus      *status;
244     rvec              mu_tot;
245     t_vcm            *vcm;
246     matrix            pcoupl_mu, M;
247     t_trxframe        rerun_fr;
248     gmx_repl_ex_t     repl_ex = NULL;
249     int               nchkpt  = 1;
250     gmx_localtop_t   *top;
251     t_mdebin         *mdebin   = NULL;
252     t_state          *state    = NULL;
253     gmx_enerdata_t   *enerd;
254     rvec             *f = NULL;
255     gmx_global_stat_t gstat;
256     gmx_update_t     *upd   = NULL;
257     t_graph          *graph = NULL;
258     gmx_groups_t     *groups;
259     gmx_ekindata_t   *ekind;
260     gmx_shellfc_t    *shellfc;
261     gmx_bool          bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
262     gmx_bool          bResetCountersHalfMaxH = FALSE;
263     gmx_bool          bTemp, bPres, bTrotter;
264     real              dvdl_constr;
265     rvec             *cbuf        = NULL;
266     int               cbuf_nalloc = 0;
267     matrix            lastbox;
268     int               lamnew  = 0;
269     /* for FEP */
270     int               nstfep = 0;
271     double            cycles;
272     real              saved_conserved_quantity = 0;
273     real              last_ekin                = 0;
274     t_extmass         MassQ;
275     int             **trotter_seq;
276     char              sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
277     int               handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
278
279
280     /* PME load balancing data for GPU kernels */
281     pme_load_balancing_t *pme_loadbal      = NULL;
282     gmx_bool              bPMETune         = FALSE;
283     gmx_bool              bPMETunePrinting = FALSE;
284
285     /* Interactive MD */
286     gmx_bool          bIMDstep = FALSE;
287
288 #ifdef GMX_FAHCORE
289     /* Temporary addition for FAHCORE checkpointing */
290     int chkpt_ret;
291 #endif
292     /* Domain decomposition could incorrectly miss a bonded
293        interaction, but checking for that requires a global
294        communication stage, which does not otherwise happen in DD
295        code. So we do that alongside the first global energy reduction
296        after a new DD is made. These variables handle whether the
297        check happens, and the result it returns. */
298     bool              shouldCheckNumberOfBondedInteractions = false;
299     int               totalNumberOfBondedInteractions       = -1;
300
301     SimulationSignals signals;
302     // Most global communnication stages don't propagate mdrun
303     // signals, and will use this object to achieve that.
304     SimulationSignaller nullSignaller(nullptr, nullptr, false, false);
305
306     /* Check for special mdrun options */
307     bRerunMD = (Flags & MD_RERUN);
308     if (Flags & MD_RESETCOUNTERSHALFWAY)
309     {
310         if (ir->nsteps > 0)
311         {
312             /* Signal to reset the counters half the simulation steps. */
313             wcycle_set_reset_counters(wcycle, ir->nsteps/2);
314         }
315         /* Signal to reset the counters halfway the simulation time. */
316         bResetCountersHalfMaxH = (max_hours > 0);
317     }
318
319     /* md-vv uses averaged full step velocities for T-control
320        md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
321        md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
322     bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
323
324     if (bRerunMD)
325     {
326         /* Since we don't know if the frames read are related in any way,
327          * rebuild the neighborlist at every step.
328          */
329         ir->nstlist       = 1;
330         ir->nstcalcenergy = 1;
331         nstglobalcomm     = 1;
332     }
333
334     nstglobalcomm   = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
335     bGStatEveryStep = (nstglobalcomm == 1);
336
337     if (bRerunMD)
338     {
339         ir->nstxout_compressed = 0;
340     }
341     groups = &top_global->groups;
342
343     if (ir->eSwapCoords != eswapNO)
344     {
345         /* Initialize ion swapping code */
346         init_swapcoords(fplog, bVerbose, ir, opt2fn_master("-swap", nfile, fnm, cr),
347                         top_global, state_global->x, state_global->box, &state_global->swapstate, cr, oenv, Flags);
348     }
349
350     /* Initial values */
351     init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
352             &(state_global->fep_state), lam0,
353             nrnb, top_global, &upd,
354             nfile, fnm, &outf, &mdebin,
355             force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags, wcycle);
356
357     clear_mat(total_vir);
358     clear_mat(pres);
359     /* Energy terms and groups */
360     snew(enerd, 1);
361     init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
362                   enerd);
363     if (DOMAINDECOMP(cr))
364     {
365         f = NULL;
366     }
367     else
368     {
369         snew(f, top_global->natoms);
370     }
371
372     /* Kinetic energy data */
373     snew(ekind, 1);
374     init_ekindata(fplog, top_global, &(ir->opts), ekind);
375     /* Copy the cos acceleration to the groups struct */
376     ekind->cosacc.cos_accel = ir->cos_accel;
377
378     gstat = global_stat_init(ir);
379
380     /* Check for polarizable models and flexible constraints */
381     shellfc = init_shell_flexcon(fplog,
382                                  top_global, n_flexible_constraints(constr),
383                                  ir->nstcalcenergy, DOMAINDECOMP(cr));
384
385     if (shellfc && ir->nstcalcenergy != 1)
386     {
387         gmx_fatal(FARGS, "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is not supported in combinations with shell particles.\nPlease make a new tpr file.", ir->nstcalcenergy);
388     }
389     if (shellfc && DOMAINDECOMP(cr))
390     {
391         gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
392     }
393
394     if (inputrecDeform(ir))
395     {
396         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
397         set_deform_reference_box(upd,
398                                  deform_init_init_step_tpx,
399                                  deform_init_box_tpx);
400         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
401     }
402
403     {
404         double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
405         if ((io > 2000) && MASTER(cr))
406         {
407             fprintf(stderr,
408                     "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
409                     io);
410         }
411     }
412
413     if (DOMAINDECOMP(cr))
414     {
415         top = dd_init_local_top(top_global);
416
417         snew(state, 1);
418         dd_init_local_state(cr->dd, state_global, state);
419     }
420     else
421     {
422         top = gmx_mtop_generate_local_top(top_global, ir->efep != efepNO);
423
424         state    = serial_init_local_state(state_global);
425
426         atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
427
428         if (vsite)
429         {
430             set_vsite_top(vsite, top, mdatoms, cr);
431         }
432
433         if (ir->ePBC != epbcNONE && !fr->bMolPBC)
434         {
435             graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
436         }
437
438         if (shellfc)
439         {
440             make_local_shells(cr, mdatoms, shellfc);
441         }
442
443         setup_bonded_threading(fr, &top->idef);
444
445         update_realloc(upd, state->nalloc);
446     }
447
448     /* Set up interactive MD (IMD) */
449     init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
450              nfile, fnm, oenv, imdport, Flags);
451
452     if (DOMAINDECOMP(cr))
453     {
454         /* Distribute the charge groups over the nodes from the master node */
455         dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
456                             state_global, top_global, ir,
457                             state, &f, mdatoms, top, fr,
458                             vsite, constr,
459                             nrnb, NULL, FALSE);
460         shouldCheckNumberOfBondedInteractions = true;
461         update_realloc(upd, state->nalloc);
462     }
463
464     update_mdatoms(mdatoms, state->lambda[efptMASS]);
465
466     startingFromCheckpoint = Flags & MD_STARTFROMCPT;
467
468     if (ir->bExpanded)
469     {
470         init_expanded_ensemble(startingFromCheckpoint, ir, &state->dfhist);
471     }
472
473     if (MASTER(cr))
474     {
475         if (startingFromCheckpoint)
476         {
477             /* Update mdebin with energy history if appending to output files */
478             if (Flags & MD_APPENDFILES)
479             {
480                 restore_energyhistory_from_state(mdebin, state_global->enerhist);
481             }
482             else
483             {
484                 /* We might have read an energy history from checkpoint,
485                  * free the allocated memory and reset the counts.
486                  */
487                 done_energyhistory(state_global->enerhist);
488                 init_energyhistory(state_global->enerhist);
489             }
490         }
491         /* Set the initial energy history in state by updating once */
492         update_energyhistory(state_global->enerhist, mdebin);
493     }
494
495     /* Initialize constraints */
496     if (constr && !DOMAINDECOMP(cr))
497     {
498         set_constraints(constr, top, ir, mdatoms, cr);
499     }
500
501     if (repl_ex_nst > 0 && MASTER(cr))
502     {
503         repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
504                                         repl_ex_nst, repl_ex_nex, repl_ex_seed);
505     }
506
507     /* PME tuning is only supported with PME for Coulomb. Is is not supported
508      * with only LJ PME, or for reruns.
509      */
510     bPMETune = ((Flags & MD_TUNEPME) && EEL_PME(fr->eeltype) && !bRerunMD &&
511                 !(Flags & MD_REPRODUCIBLE));
512     if (bPMETune)
513     {
514         pme_loadbal_init(&pme_loadbal, cr, fplog, ir, state->box,
515                          fr->ic, fr->pmedata, use_GPU(fr->nbv),
516                          &bPMETunePrinting);
517     }
518
519     if (!ir->bContinuation && !bRerunMD)
520     {
521         if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
522         {
523             /* Set the velocities of frozen particles to zero */
524             for (i = 0; i < mdatoms->homenr; i++)
525             {
526                 for (m = 0; m < DIM; m++)
527                 {
528                     if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
529                     {
530                         state->v[i][m] = 0;
531                     }
532                 }
533             }
534         }
535
536         if (constr)
537         {
538             /* Constrain the initial coordinates and velocities */
539             do_constrain_first(fplog, constr, ir, mdatoms, state,
540                                cr, nrnb, fr, top);
541         }
542         if (vsite)
543         {
544             /* Construct the virtual sites for the initial configuration */
545             construct_vsites(vsite, state->x, ir->delta_t, NULL,
546                              top->idef.iparams, top->idef.il,
547                              fr->ePBC, fr->bMolPBC, cr, state->box);
548         }
549     }
550
551     if (ir->efep != efepNO)
552     {
553         /* Set free energy calculation frequency as the greatest common
554          * denominator of nstdhdl and repl_ex_nst. */
555         nstfep = ir->fepvals->nstdhdl;
556         if (ir->bExpanded)
557         {
558             nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
559         }
560         if (repl_ex_nst > 0)
561         {
562             nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
563         }
564     }
565
566     /* Be REALLY careful about what flags you set here. You CANNOT assume
567      * this is the first step, since we might be restarting from a checkpoint,
568      * and in that case we should not do any modifications to the state.
569      */
570     bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
571
572     if (Flags & MD_READ_EKIN)
573     {
574         restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
575     }
576
577     cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
578                   | (bStopCM ? CGLO_STOPCM : 0)
579                   | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
580                   | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
581                   | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
582
583     bSumEkinhOld = FALSE;
584     compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
585                     NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
586                     constr, &nullSignaller, state->box,
587                     &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags
588                     | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
589     checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
590                                     top_global, top, state,
591                                     &shouldCheckNumberOfBondedInteractions);
592     if (ir->eI == eiVVAK)
593     {
594         /* a second call to get the half step temperature initialized as well */
595         /* we do the same call as above, but turn the pressure off -- internally to
596            compute_globals, this is recognized as a velocity verlet half-step
597            kinetic energy calculation.  This minimized excess variables, but
598            perhaps loses some logic?*/
599
600         compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
601                         NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
602                         constr, &nullSignaller, state->box,
603                         NULL, &bSumEkinhOld,
604                         cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
605     }
606
607     /* Calculate the initial half step temperature, and save the ekinh_old */
608     if (!(Flags & MD_STARTFROMCPT))
609     {
610         for (i = 0; (i < ir->opts.ngtc); i++)
611         {
612             copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
613         }
614     }
615     if (ir->eI != eiVV)
616     {
617         enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
618                                      and there is no previous step */
619     }
620
621     /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
622        temperature control */
623     trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
624
625     if (MASTER(cr))
626     {
627         if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
628         {
629             fprintf(fplog,
630                     "RMS relative constraint deviation after constraining: %.2e\n",
631                     constr_rmsd(constr));
632         }
633         if (EI_STATE_VELOCITY(ir->eI))
634         {
635             fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
636         }
637         if (bRerunMD)
638         {
639             fprintf(stderr, "starting md rerun '%s', reading coordinates from"
640                     " input trajectory '%s'\n\n",
641                     *(top_global->name), opt2fn("-rerun", nfile, fnm));
642             if (bVerbose)
643             {
644                 fprintf(stderr, "Calculated time to finish depends on nsteps from "
645                         "run input file,\nwhich may not correspond to the time "
646                         "needed to process input trajectory.\n\n");
647             }
648         }
649         else
650         {
651             char tbuf[20];
652             fprintf(stderr, "starting mdrun '%s'\n",
653                     *(top_global->name));
654             if (ir->nsteps >= 0)
655             {
656                 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
657             }
658             else
659             {
660                 sprintf(tbuf, "%s", "infinite");
661             }
662             if (ir->init_step > 0)
663             {
664                 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
665                         gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
666                         gmx_step_str(ir->init_step, sbuf2),
667                         ir->init_step*ir->delta_t);
668             }
669             else
670             {
671                 fprintf(stderr, "%s steps, %s ps.\n",
672                         gmx_step_str(ir->nsteps, sbuf), tbuf);
673             }
674         }
675         fprintf(fplog, "\n");
676     }
677
678     walltime_accounting_start(walltime_accounting);
679     wallcycle_start(wcycle, ewcRUN);
680     print_start(fplog, cr, walltime_accounting, "mdrun");
681
682     /* safest point to do file checkpointing is here.  More general point would be immediately before integrator call */
683 #ifdef GMX_FAHCORE
684     chkpt_ret = fcCheckPointParallel( cr->nodeid,
685                                       NULL, 0);
686     if (chkpt_ret == 0)
687     {
688         gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
689     }
690 #endif
691
692     /***********************************************************
693      *
694      *             Loop over MD steps
695      *
696      ************************************************************/
697
698     /* if rerunMD then read coordinates and velocities from input trajectory */
699     if (bRerunMD)
700     {
701         if (getenv("GMX_FORCE_UPDATE"))
702         {
703             bForceUpdate = TRUE;
704         }
705
706         rerun_fr.natoms = 0;
707         if (MASTER(cr))
708         {
709             bLastStep = !read_first_frame(oenv, &status,
710                                           opt2fn("-rerun", nfile, fnm),
711                                           &rerun_fr, TRX_NEED_X | TRX_READ_V);
712             if (rerun_fr.natoms != top_global->natoms)
713             {
714                 gmx_fatal(FARGS,
715                           "Number of atoms in trajectory (%d) does not match the "
716                           "run input file (%d)\n",
717                           rerun_fr.natoms, top_global->natoms);
718             }
719             if (ir->ePBC != epbcNONE)
720             {
721                 if (!rerun_fr.bBox)
722                 {
723                     gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f does not contain a box, while pbc is used", rerun_fr.step, rerun_fr.time);
724                 }
725                 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
726                 {
727                     gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
728                 }
729             }
730         }
731
732         if (PAR(cr))
733         {
734             rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
735         }
736
737         if (ir->ePBC != epbcNONE)
738         {
739             /* Set the shift vectors.
740              * Necessary here when have a static box different from the tpr box.
741              */
742             calc_shifts(rerun_fr.box, fr->shift_vec);
743         }
744     }
745
746     /* loop over MD steps or if rerunMD to end of input trajectory */
747     bFirstStep = TRUE;
748     /* Skip the first Nose-Hoover integration when we get the state from tpx */
749     bInitStep        = !startingFromCheckpoint || EI_VV(ir->eI);
750     bSumEkinhOld     = FALSE;
751     bExchanged       = FALSE;
752     bNeedRepartition = FALSE;
753     // TODO This implementation of ensemble orientation restraints is nasty because
754     // a user can't just do multi-sim with single-sim orientation restraints.
755     bUsingEnsembleRestraints = (fcd->disres.nsystems > 1) || (cr->ms && fcd->orires.nr);
756
757     {
758         // Replica exchange and ensemble restraints need all
759         // simulations to remain synchronized, so they need
760         // checkpoints and stop conditions to act on the same step, so
761         // the propagation of such signals must take place between
762         // simulations, not just within simulations.
763         bool checkpointIsLocal    = (repl_ex_nst <= 0) && !bUsingEnsembleRestraints;
764         bool stopConditionIsLocal = (repl_ex_nst <= 0) && !bUsingEnsembleRestraints;
765         bool resetCountersIsLocal = true;
766         signals[eglsCHKPT]         = SimulationSignal(checkpointIsLocal);
767         signals[eglsSTOPCOND]      = SimulationSignal(stopConditionIsLocal);
768         signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
769     }
770
771     step     = ir->init_step;
772     step_rel = 0;
773
774     // TODO extract this to new multi-simulation module
775     if (MASTER(cr) && MULTISIM(cr) && (repl_ex_nst <= 0 ))
776     {
777         if (!multisim_int_all_are_equal(cr->ms, ir->nsteps))
778         {
779             md_print_info(cr, fplog,
780                           "Note: The number of steps is not consistent across multi simulations,\n"
781                           "but we are proceeding anyway!\n");
782         }
783         if (!multisim_int_all_are_equal(cr->ms, ir->init_step))
784         {
785             md_print_info(cr, fplog,
786                           "Note: The initial step is not consistent across multi simulations,\n"
787                           "but we are proceeding anyway!\n");
788         }
789     }
790
791     /* and stop now if we should */
792     bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
793     while (!bLastStep)
794     {
795
796         /* Determine if this is a neighbor search step */
797         bNStList = (ir->nstlist > 0  && step % ir->nstlist == 0);
798
799         if (bPMETune && bNStList)
800         {
801             /* PME grid + cut-off optimization with GPUs or PME nodes */
802             pme_loadbal_do(pme_loadbal, cr,
803                            (bVerbose && MASTER(cr)) ? stderr : NULL,
804                            fplog,
805                            ir, fr, state,
806                            wcycle,
807                            step, step_rel,
808                            &bPMETunePrinting);
809         }
810
811         wallcycle_start(wcycle, ewcSTEP);
812
813         if (bRerunMD)
814         {
815             if (rerun_fr.bStep)
816             {
817                 step     = rerun_fr.step;
818                 step_rel = step - ir->init_step;
819             }
820             if (rerun_fr.bTime)
821             {
822                 t = rerun_fr.time;
823             }
824             else
825             {
826                 t = step;
827             }
828         }
829         else
830         {
831             bLastStep = (step_rel == ir->nsteps);
832             t         = t0 + step*ir->delta_t;
833         }
834
835         // TODO Refactor this, so that nstfep does not need a default value of zero
836         if (ir->efep != efepNO || ir->bSimTemp)
837         {
838             /* find and set the current lambdas.  If rerunning, we either read in a state, or a lambda value,
839                requiring different logic. */
840
841             set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
842             bDoDHDL      = do_per_step(step, ir->fepvals->nstdhdl);
843             bDoFEP       = ((ir->efep != efepNO) && do_per_step(step, nstfep));
844             bDoExpanded  = (do_per_step(step, ir->expandedvals->nstexpanded)
845                             && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
846         }
847
848         bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
849                      do_per_step(step, repl_ex_nst));
850
851         if (bSimAnn)
852         {
853             update_annealing_target_temp(ir, t, upd);
854         }
855
856         if (bRerunMD)
857         {
858             if (!DOMAINDECOMP(cr) || MASTER(cr))
859             {
860                 for (i = 0; i < state_global->natoms; i++)
861                 {
862                     copy_rvec(rerun_fr.x[i], state_global->x[i]);
863                 }
864                 if (rerun_fr.bV)
865                 {
866                     for (i = 0; i < state_global->natoms; i++)
867                     {
868                         copy_rvec(rerun_fr.v[i], state_global->v[i]);
869                     }
870                 }
871                 else
872                 {
873                     for (i = 0; i < state_global->natoms; i++)
874                     {
875                         clear_rvec(state_global->v[i]);
876                     }
877                     if (bRerunWarnNoV)
878                     {
879                         fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
880                                 "         Ekin, temperature and pressure are incorrect,\n"
881                                 "         the virial will be incorrect when constraints are present.\n"
882                                 "\n");
883                         bRerunWarnNoV = FALSE;
884                     }
885                 }
886             }
887             copy_mat(rerun_fr.box, state_global->box);
888             copy_mat(state_global->box, state->box);
889
890             if (vsite && (Flags & MD_RERUN_VSITE))
891             {
892                 if (DOMAINDECOMP(cr))
893                 {
894                     gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
895                 }
896                 if (graph)
897                 {
898                     /* Following is necessary because the graph may get out of sync
899                      * with the coordinates if we only have every N'th coordinate set
900                      */
901                     mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
902                     shift_self(graph, state->box, state->x);
903                 }
904                 construct_vsites(vsite, state->x, ir->delta_t, state->v,
905                                  top->idef.iparams, top->idef.il,
906                                  fr->ePBC, fr->bMolPBC, cr, state->box);
907                 if (graph)
908                 {
909                     unshift_self(graph, state->box, state->x);
910                 }
911             }
912         }
913
914         /* Stop Center of Mass motion */
915         bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
916
917         if (bRerunMD)
918         {
919             /* for rerun MD always do Neighbour Searching */
920             bNS      = (bFirstStep || ir->nstlist != 0);
921         }
922         else
923         {
924             /* Determine whether or not to do Neighbour Searching */
925             bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
926         }
927
928         /* < 0 means stop at next step, > 0 means stop at next NS step */
929         if ( (signals[eglsSTOPCOND].set < 0) ||
930              ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
931         {
932             bLastStep = TRUE;
933         }
934
935         /* Determine whether or not to update the Born radii if doing GB */
936         bBornRadii = bFirstStep;
937         if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
938         {
939             bBornRadii = TRUE;
940         }
941
942         /* do_log triggers energy and virial calculation. Because this leads
943          * to different code paths, forces can be different. Thus for exact
944          * continuation we should avoid extra log output.
945          * Note that the || bLastStep can result in non-exact continuation
946          * beyond the last step. But we don't consider that to be an issue.
947          */
948         do_log     = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
949         do_verbose = bVerbose &&
950             (step % stepout == 0 || bFirstStep || bLastStep || bRerunMD);
951
952         if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
953         {
954             if (bRerunMD)
955             {
956                 bMasterState = TRUE;
957             }
958             else
959             {
960                 bMasterState = FALSE;
961                 /* Correct the new box if it is too skewed */
962                 if (inputrecDynamicBox(ir))
963                 {
964                     if (correct_box(fplog, step, state->box, graph))
965                     {
966                         bMasterState = TRUE;
967                     }
968                 }
969                 if (DOMAINDECOMP(cr) && bMasterState)
970                 {
971                     dd_collect_state(cr->dd, state, state_global);
972                 }
973             }
974
975             if (DOMAINDECOMP(cr))
976             {
977                 /* Repartition the domain decomposition */
978                 dd_partition_system(fplog, step, cr,
979                                     bMasterState, nstglobalcomm,
980                                     state_global, top_global, ir,
981                                     state, &f, mdatoms, top, fr,
982                                     vsite, constr,
983                                     nrnb, wcycle,
984                                     do_verbose && !bPMETunePrinting);
985                 shouldCheckNumberOfBondedInteractions = true;
986                 update_realloc(upd, state->nalloc);
987             }
988         }
989
990         if (MASTER(cr) && do_log)
991         {
992             print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
993         }
994
995         if (ir->efep != efepNO)
996         {
997             update_mdatoms(mdatoms, state->lambda[efptMASS]);
998         }
999
1000         if ((bRerunMD && rerun_fr.bV) || bExchanged)
1001         {
1002
1003             /* We need the kinetic energy at minus the half step for determining
1004              * the full step kinetic energy and possibly for T-coupling.*/
1005             /* This may not be quite working correctly yet . . . . */
1006             compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1007                             wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1008                             constr, &nullSignaller, state->box,
1009                             &totalNumberOfBondedInteractions, &bSumEkinhOld,
1010                             CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1011             checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1012                                             top_global, top, state,
1013                                             &shouldCheckNumberOfBondedInteractions);
1014         }
1015         clear_mat(force_vir);
1016
1017         /* We write a checkpoint at this MD step when:
1018          * either at an NS step when we signalled through gs,
1019          * or at the last step (but not when we do not want confout),
1020          * but never at the first step or with rerun.
1021          */
1022         bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1023                  (bLastStep && (Flags & MD_CONFOUT))) &&
1024                 step > ir->init_step && !bRerunMD);
1025         if (bCPT)
1026         {
1027             signals[eglsCHKPT].set = 0;
1028         }
1029
1030         /* Determine the energy and pressure:
1031          * at nstcalcenergy steps and at energy output steps (set below).
1032          */
1033         if (EI_VV(ir->eI) && (!bInitStep))
1034         {
1035             /* for vv, the first half of the integration actually corresponds
1036                to the previous step.  bCalcEner is only required to be evaluated on the 'next' step,
1037                but the virial needs to be calculated on both the current step and the 'next' step. Future
1038                reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1039
1040             /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1041             bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1042             bCalcVir      = bCalcEnerStep ||
1043                 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1044         }
1045         else
1046         {
1047             bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1048             bCalcVir      = bCalcEnerStep ||
1049                 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1050         }
1051         bCalcEner = bCalcEnerStep;
1052
1053         do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1054
1055         if (do_ene || do_log || bDoReplEx)
1056         {
1057             bCalcVir  = TRUE;
1058             bCalcEner = TRUE;
1059         }
1060
1061         /* Do we need global communication ? */
1062         bGStat = (bCalcVir || bCalcEner || bStopCM ||
1063                   do_per_step(step, nstglobalcomm) ||
1064                   (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1065
1066         force_flags = (GMX_FORCE_STATECHANGED |
1067                        ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1068                        GMX_FORCE_ALLFORCES |
1069                        (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1070                        (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1071                        (bDoFEP ? GMX_FORCE_DHDL : 0)
1072                        );
1073
1074         if (shellfc)
1075         {
1076             /* Now is the time to relax the shells */
1077             relax_shell_flexcon(fplog, cr, bVerbose, step,
1078                                 ir, bNS, force_flags, top,
1079                                 constr, enerd, fcd,
1080                                 state, f, force_vir, mdatoms,
1081                                 nrnb, wcycle, graph, groups,
1082                                 shellfc, fr, bBornRadii, t, mu_tot,
1083                                 vsite, mdoutf_get_fp_field(outf));
1084         }
1085         else
1086         {
1087             /* The coordinates (x) are shifted (to get whole molecules)
1088              * in do_force.
1089              * This is parallellized as well, and does communication too.
1090              * Check comments in sim_util.c
1091              */
1092             do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1093                      state->box, state->x, &state->hist,
1094                      f, force_vir, mdatoms, enerd, fcd,
1095                      state->lambda, graph,
1096                      fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1097                      (bNS ? GMX_FORCE_NS : 0) | force_flags);
1098         }
1099
1100         if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1101         /*  ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1102         {
1103             rvec *vbuf = NULL;
1104
1105             wallcycle_start(wcycle, ewcUPDATE);
1106             if (ir->eI == eiVV && bInitStep)
1107             {
1108                 /* if using velocity verlet with full time step Ekin,
1109                  * take the first half step only to compute the
1110                  * virial for the first step. From there,
1111                  * revert back to the initial coordinates
1112                  * so that the input is actually the initial step.
1113                  */
1114                 snew(vbuf, state->natoms);
1115                 copy_rvecn(state->v, vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1116             }
1117             else
1118             {
1119                 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1120                 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1121             }
1122
1123             update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1124                           ekind, M, upd, etrtVELOCITY1,
1125                           cr, constr);
1126
1127             if (!bRerunMD || rerun_fr.bV || bForceUpdate)         /* Why is rerun_fr.bV here?  Unclear. */
1128             {
1129                 wallcycle_stop(wcycle, ewcUPDATE);
1130                 update_constraints(fplog, step, NULL, ir, mdatoms,
1131                                    state, fr->bMolPBC, graph, f,
1132                                    &top->idef, shake_vir,
1133                                    cr, nrnb, wcycle, upd, constr,
1134                                    TRUE, bCalcVir);
1135                 wallcycle_start(wcycle, ewcUPDATE);
1136             }
1137             else if (graph)
1138             {
1139                 /* Need to unshift here if a do_force has been
1140                    called in the previous step */
1141                 unshift_self(graph, state->box, state->x);
1142             }
1143             /* if VV, compute the pressure and constraints */
1144             /* For VV2, we strictly only need this if using pressure
1145              * control, but we really would like to have accurate pressures
1146              * printed out.
1147              * Think about ways around this in the future?
1148              * For now, keep this choice in comments.
1149              */
1150             /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1151             /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1152             bPres = TRUE;
1153             bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1154             if (bCalcEner && ir->eI == eiVVAK)
1155             {
1156                 bSumEkinhOld = TRUE;
1157             }
1158             /* for vv, the first half of the integration actually corresponds to the previous step.
1159                So we need information from the last step in the first half of the integration */
1160             if (bGStat || do_per_step(step-1, nstglobalcomm))
1161             {
1162                 wallcycle_stop(wcycle, ewcUPDATE);
1163                 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1164                                 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1165                                 constr, &nullSignaller, state->box,
1166                                 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1167                                 (bGStat ? CGLO_GSTAT : 0)
1168                                 | CGLO_ENERGY
1169                                 | (bTemp ? CGLO_TEMPERATURE : 0)
1170                                 | (bPres ? CGLO_PRESSURE : 0)
1171                                 | (bPres ? CGLO_CONSTRAINT : 0)
1172                                 | (bStopCM ? CGLO_STOPCM : 0)
1173                                 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1174                                 | CGLO_SCALEEKIN
1175                                 );
1176                 /* explanation of above:
1177                    a) We compute Ekin at the full time step
1178                    if 1) we are using the AveVel Ekin, and it's not the
1179                    initial step, or 2) if we are using AveEkin, but need the full
1180                    time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1181                    b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1182                    EkinAveVel because it's needed for the pressure */
1183                 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1184                                                 top_global, top, state,
1185                                                 &shouldCheckNumberOfBondedInteractions);
1186                 wallcycle_start(wcycle, ewcUPDATE);
1187             }
1188             /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1189             if (!bInitStep)
1190             {
1191                 if (bTrotter)
1192                 {
1193                     m_add(force_vir, shake_vir, total_vir);     /* we need the un-dispersion corrected total vir here */
1194                     trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1195
1196                     copy_mat(shake_vir, state->svir_prev);
1197                     copy_mat(force_vir, state->fvir_prev);
1198                     if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1199                     {
1200                         /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1201                         enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1202                         enerd->term[F_EKIN] = trace(ekind->ekin);
1203                     }
1204                 }
1205                 else if (bExchanged)
1206                 {
1207                     wallcycle_stop(wcycle, ewcUPDATE);
1208                     /* We need the kinetic energy at minus the half step for determining
1209                      * the full step kinetic energy and possibly for T-coupling.*/
1210                     /* This may not be quite working correctly yet . . . . */
1211                     compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1212                                     wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1213                                     constr, &nullSignaller, state->box,
1214                                     NULL, &bSumEkinhOld,
1215                                     CGLO_GSTAT | CGLO_TEMPERATURE);
1216                     wallcycle_start(wcycle, ewcUPDATE);
1217                 }
1218             }
1219             /* if it's the initial step, we performed this first step just to get the constraint virial */
1220             if (ir->eI == eiVV && bInitStep)
1221             {
1222                 copy_rvecn(vbuf, state->v, 0, state->natoms);
1223                 sfree(vbuf);
1224             }
1225             wallcycle_stop(wcycle, ewcUPDATE);
1226         }
1227
1228         /* compute the conserved quantity */
1229         if (EI_VV(ir->eI))
1230         {
1231             saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1232             if (ir->eI == eiVV)
1233             {
1234                 last_ekin = enerd->term[F_EKIN];
1235             }
1236             if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1237             {
1238                 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1239             }
1240             /* sum up the foreign energy and dhdl terms for vv.  currently done every step so that dhdl is correct in the .edr */
1241             if (ir->efep != efepNO && !bRerunMD)
1242             {
1243                 sum_dhdl(enerd, state->lambda, ir->fepvals);
1244             }
1245         }
1246
1247         /* ########  END FIRST UPDATE STEP  ############## */
1248         /* ########  If doing VV, we now have v(dt) ###### */
1249         if (bDoExpanded)
1250         {
1251             /* perform extended ensemble sampling in lambda - we don't
1252                actually move to the new state before outputting
1253                statistics, but if performing simulated tempering, we
1254                do update the velocities and the tau_t. */
1255
1256             lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1257             /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1258             copy_df_history(&state_global->dfhist, &state->dfhist);
1259         }
1260
1261         /* Now we have the energies and forces corresponding to the
1262          * coordinates at time t. We must output all of this before
1263          * the update.
1264          */
1265         do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1266                                  ir, state, state_global, top_global, fr,
1267                                  outf, mdebin, ekind, f,
1268                                  &nchkpt,
1269                                  bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1270                                  bSumEkinhOld);
1271         /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1272         bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1273
1274         /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1275         if (startingFromCheckpoint && bTrotter)
1276         {
1277             copy_mat(state->svir_prev, shake_vir);
1278             copy_mat(state->fvir_prev, force_vir);
1279         }
1280
1281         elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1282
1283         /* Check whether everything is still allright */
1284         if (((int)gmx_get_stop_condition() > handled_stop_condition)
1285 #if GMX_THREAD_MPI
1286             && MASTER(cr)
1287 #endif
1288             )
1289         {
1290             int nsteps_stop = -1;
1291
1292             /* this just makes signals[].sig compatible with the hack
1293                of sending signals around by MPI_Reduce together with
1294                other floats */
1295             if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1296             {
1297                 signals[eglsSTOPCOND].sig = 1;
1298                 nsteps_stop               = std::max(ir->nstlist, 2*nstglobalcomm);
1299             }
1300             if (gmx_get_stop_condition() == gmx_stop_cond_next)
1301             {
1302                 signals[eglsSTOPCOND].sig = -1;
1303                 nsteps_stop               = nstglobalcomm + 1;
1304             }
1305             if (fplog)
1306             {
1307                 fprintf(fplog,
1308                         "\n\nReceived the %s signal, stopping within %d steps\n\n",
1309                         gmx_get_signal_name(), nsteps_stop);
1310                 fflush(fplog);
1311             }
1312             fprintf(stderr,
1313                     "\n\nReceived the %s signal, stopping within %d steps\n\n",
1314                     gmx_get_signal_name(), nsteps_stop);
1315             fflush(stderr);
1316             handled_stop_condition = (int)gmx_get_stop_condition();
1317         }
1318         else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1319                  (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1320                  signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1321         {
1322             /* Signal to terminate the run */
1323             signals[eglsSTOPCOND].sig = 1;
1324             if (fplog)
1325             {
1326                 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1327             }
1328             fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1329         }
1330
1331         if (bResetCountersHalfMaxH && MASTER(cr) &&
1332             elapsed_time > max_hours*60.0*60.0*0.495)
1333         {
1334             /* Set flag that will communicate the signal to all ranks in the simulation */
1335             signals[eglsRESETCOUNTERS].sig = 1;
1336         }
1337
1338         /* In parallel we only have to check for checkpointing in steps
1339          * where we do global communication,
1340          *  otherwise the other nodes don't know.
1341          */
1342         if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1343                            cpt_period >= 0 &&
1344                            (cpt_period == 0 ||
1345                             elapsed_time >= nchkpt*cpt_period*60.0)) &&
1346             signals[eglsCHKPT].set == 0)
1347         {
1348             signals[eglsCHKPT].sig = 1;
1349         }
1350
1351         /* #########   START SECOND UPDATE STEP ################# */
1352
1353         /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1354            in preprocessing */
1355
1356         if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1357         {
1358             gmx_bool bIfRandomize;
1359             bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1360             /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1361             if (constr && bIfRandomize)
1362             {
1363                 update_constraints(fplog, step, NULL, ir, mdatoms,
1364                                    state, fr->bMolPBC, graph, f,
1365                                    &top->idef, tmp_vir,
1366                                    cr, nrnb, wcycle, upd, constr,
1367                                    TRUE, bCalcVir);
1368             }
1369         }
1370         /* Box is changed in update() when we do pressure coupling,
1371          * but we should still use the old box for energy corrections and when
1372          * writing it to the energy file, so it matches the trajectory files for
1373          * the same timestep above. Make a copy in a separate array.
1374          */
1375         copy_mat(state->box, lastbox);
1376
1377         dvdl_constr = 0;
1378
1379         if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1380         {
1381             wallcycle_start(wcycle, ewcUPDATE);
1382             /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1383             if (bTrotter)
1384             {
1385                 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1386                 /* We can only do Berendsen coupling after we have summed
1387                  * the kinetic energy or virial. Since the happens
1388                  * in global_state after update, we should only do it at
1389                  * step % nstlist = 1 with bGStatEveryStep=FALSE.
1390                  */
1391             }
1392             else
1393             {
1394                 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1395                 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1396             }
1397
1398             if (EI_VV(ir->eI))
1399             {
1400                 /* velocity half-step update */
1401                 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1402                               ekind, M, upd, etrtVELOCITY2,
1403                               cr, constr);
1404             }
1405
1406             /* Above, initialize just copies ekinh into ekin,
1407              * it doesn't copy position (for VV),
1408              * and entire integrator for MD.
1409              */
1410
1411             if (ir->eI == eiVVAK)
1412             {
1413                 /* We probably only need md->homenr, not state->natoms */
1414                 if (state->natoms > cbuf_nalloc)
1415                 {
1416                     cbuf_nalloc = state->natoms;
1417                     srenew(cbuf, cbuf_nalloc);
1418                 }
1419                 copy_rvecn(state->x, cbuf, 0, state->natoms);
1420             }
1421
1422             update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1423                           ekind, M, upd, etrtPOSITION, cr, constr);
1424             wallcycle_stop(wcycle, ewcUPDATE);
1425
1426             update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1427                                fr->bMolPBC, graph, f,
1428                                &top->idef, shake_vir,
1429                                cr, nrnb, wcycle, upd, constr,
1430                                FALSE, bCalcVir);
1431
1432             if (ir->eI == eiVVAK)
1433             {
1434                 /* erase F_EKIN and F_TEMP here? */
1435                 /* just compute the kinetic energy at the half step to perform a trotter step */
1436                 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1437                                 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1438                                 constr, &nullSignaller, lastbox,
1439                                 NULL, &bSumEkinhOld,
1440                                 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1441                                 );
1442                 wallcycle_start(wcycle, ewcUPDATE);
1443                 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1444                 /* now we know the scaling, we can compute the positions again again */
1445                 copy_rvecn(cbuf, state->x, 0, state->natoms);
1446
1447                 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1448                               ekind, M, upd, etrtPOSITION, cr, constr);
1449                 wallcycle_stop(wcycle, ewcUPDATE);
1450
1451                 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1452                 /* are the small terms in the shake_vir here due
1453                  * to numerical errors, or are they important
1454                  * physically? I'm thinking they are just errors, but not completely sure.
1455                  * For now, will call without actually constraining, constr=NULL*/
1456                 update_constraints(fplog, step, NULL, ir, mdatoms,
1457                                    state, fr->bMolPBC, graph, f,
1458                                    &top->idef, tmp_vir,
1459                                    cr, nrnb, wcycle, upd, NULL,
1460                                    FALSE, bCalcVir);
1461             }
1462             if (EI_VV(ir->eI))
1463             {
1464                 /* this factor or 2 correction is necessary
1465                    because half of the constraint force is removed
1466                    in the vv step, so we have to double it.  See
1467                    the Redmine issue #1255.  It is not yet clear
1468                    if the factor of 2 is exact, or just a very
1469                    good approximation, and this will be
1470                    investigated.  The next step is to see if this
1471                    can be done adding a dhdl contribution from the
1472                    rattle step, but this is somewhat more
1473                    complicated with the current code. Will be
1474                    investigated, hopefully for 4.6.3. However,
1475                    this current solution is much better than
1476                    having it completely wrong.
1477                  */
1478                 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1479             }
1480             else
1481             {
1482                 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1483             }
1484         }
1485         else if (graph)
1486         {
1487             /* Need to unshift here */
1488             unshift_self(graph, state->box, state->x);
1489         }
1490
1491         if (vsite != NULL)
1492         {
1493             wallcycle_start(wcycle, ewcVSITECONSTR);
1494             if (graph != NULL)
1495             {
1496                 shift_self(graph, state->box, state->x);
1497             }
1498             construct_vsites(vsite, state->x, ir->delta_t, state->v,
1499                              top->idef.iparams, top->idef.il,
1500                              fr->ePBC, fr->bMolPBC, cr, state->box);
1501
1502             if (graph != NULL)
1503             {
1504                 unshift_self(graph, state->box, state->x);
1505             }
1506             wallcycle_stop(wcycle, ewcVSITECONSTR);
1507         }
1508
1509         /* ############## IF NOT VV, Calculate globals HERE  ############ */
1510         /* With Leap-Frog we can skip compute_globals at
1511          * non-communication steps, but we need to calculate
1512          * the kinetic energy one step before communication.
1513          */
1514         {
1515             // Organize to do inter-simulation signalling on steps if
1516             // and when algorithms require it.
1517             bool doInterSimSignal = (!bFirstStep && bDoReplEx) || bUsingEnsembleRestraints;
1518
1519             if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1520             {
1521                 // Since we're already communicating at this step, we
1522                 // can propagate intra-simulation signals. Note that
1523                 // check_nstglobalcomm has the responsibility for
1524                 // choosing the value of nstglobalcomm that is one way
1525                 // bGStat becomes true, so we can't get into a
1526                 // situation where e.g. checkpointing can't be
1527                 // signalled.
1528                 bool                doIntraSimSignal = true;
1529                 SimulationSignaller signaller(&signals, cr, doInterSimSignal, doIntraSimSignal);
1530
1531                 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1532                                 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1533                                 constr, &signaller,
1534                                 lastbox,
1535                                 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1536                                 (bGStat ? CGLO_GSTAT : 0)
1537                                 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1538                                 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1539                                 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1540                                 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1541                                 | CGLO_CONSTRAINT
1542                                 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1543                                 );
1544                 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1545                                                 top_global, top, state,
1546                                                 &shouldCheckNumberOfBondedInteractions);
1547             }
1548         }
1549
1550         /* #############  END CALC EKIN AND PRESSURE ################# */
1551
1552         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1553            the virial that should probably be addressed eventually. state->veta has better properies,
1554            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1555            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
1556
1557         if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1558         {
1559             /* Sum up the foreign energy and dhdl terms for md and sd.
1560                Currently done every step so that dhdl is correct in the .edr */
1561             sum_dhdl(enerd, state->lambda, ir->fepvals);
1562         }
1563         update_box(fplog, step, ir, mdatoms, state, f,
1564                    pcoupl_mu, nrnb, upd);
1565
1566         /* ################# END UPDATE STEP 2 ################# */
1567         /* #### We now have r(t+dt) and v(t+dt/2)  ############# */
1568
1569         /* The coordinates (x) were unshifted in update */
1570         if (!bGStat)
1571         {
1572             /* We will not sum ekinh_old,
1573              * so signal that we still have to do it.
1574              */
1575             bSumEkinhOld = TRUE;
1576         }
1577
1578         /* #########  BEGIN PREPARING EDR OUTPUT  ###########  */
1579
1580         /* use the directly determined last velocity, not actually the averaged half steps */
1581         if (bTrotter && ir->eI == eiVV)
1582         {
1583             enerd->term[F_EKIN] = last_ekin;
1584         }
1585         enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1586
1587         if (EI_VV(ir->eI))
1588         {
1589             enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1590         }
1591         else
1592         {
1593             enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1594         }
1595         /* #########  END PREPARING EDR OUTPUT  ###########  */
1596
1597         /* Output stuff */
1598         if (MASTER(cr))
1599         {
1600             if (fplog && do_log && bDoExpanded)
1601             {
1602                 /* only needed if doing expanded ensemble */
1603                 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1604                                           &state_global->dfhist, state->fep_state, ir->nstlog, step);
1605             }
1606             if (bCalcEner)
1607             {
1608                 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1609                            t, mdatoms->tmass, enerd, state,
1610                            ir->fepvals, ir->expandedvals, lastbox,
1611                            shake_vir, force_vir, total_vir, pres,
1612                            ekind, mu_tot, constr);
1613             }
1614             else
1615             {
1616                 upd_mdebin_step(mdebin);
1617             }
1618
1619             gmx_bool do_dr  = do_per_step(step, ir->nstdisreout);
1620             gmx_bool do_or  = do_per_step(step, ir->nstorireout);
1621
1622             print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1623                        step, t,
1624                        eprNORMAL, mdebin, fcd, groups, &(ir->opts));
1625
1626             if (ir->bPull)
1627             {
1628                 pull_print_output(ir->pull_work, step, t);
1629             }
1630
1631             if (do_per_step(step, ir->nstlog))
1632             {
1633                 if (fflush(fplog) != 0)
1634                 {
1635                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1636                 }
1637             }
1638         }
1639         if (bDoExpanded)
1640         {
1641             /* Have to do this part _after_ outputting the logfile and the edr file */
1642             /* Gets written into the state at the beginning of next loop*/
1643             state->fep_state = lamnew;
1644         }
1645         /* Print the remaining wall clock time for the run */
1646         if (MULTIMASTER(cr) &&
1647             (do_verbose || gmx_got_usr_signal()) &&
1648             !bPMETunePrinting)
1649         {
1650             if (shellfc)
1651             {
1652                 fprintf(stderr, "\n");
1653             }
1654             print_time(stderr, walltime_accounting, step, ir, cr);
1655         }
1656
1657         /* Ion/water position swapping.
1658          * Not done in last step since trajectory writing happens before this call
1659          * in the MD loop and exchanges would be lost anyway. */
1660         bNeedRepartition = FALSE;
1661         if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1662             do_per_step(step, ir->swap->nstswap))
1663         {
1664             bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1665                                              bRerunMD ? rerun_fr.x   : state->x,
1666                                              bRerunMD ? rerun_fr.box : state->box,
1667                                              top_global, MASTER(cr) && bVerbose, bRerunMD);
1668
1669             if (bNeedRepartition && DOMAINDECOMP(cr))
1670             {
1671                 dd_collect_state(cr->dd, state, state_global);
1672             }
1673         }
1674
1675         /* Replica exchange */
1676         bExchanged = FALSE;
1677         if (bDoReplEx)
1678         {
1679             bExchanged = replica_exchange(fplog, cr, repl_ex,
1680                                           state_global, enerd,
1681                                           state, step, t);
1682         }
1683
1684         if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1685         {
1686             dd_partition_system(fplog, step, cr, TRUE, 1,
1687                                 state_global, top_global, ir,
1688                                 state, &f, mdatoms, top, fr,
1689                                 vsite, constr,
1690                                 nrnb, wcycle, FALSE);
1691             shouldCheckNumberOfBondedInteractions = true;
1692             update_realloc(upd, state->nalloc);
1693         }
1694
1695         bFirstStep             = FALSE;
1696         bInitStep              = FALSE;
1697         startingFromCheckpoint = FALSE;
1698
1699         /* #######  SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1700         /* With all integrators, except VV, we need to retain the pressure
1701          * at the current step for coupling at the next step.
1702          */
1703         if ((state->flags & (1<<estPRES_PREV)) &&
1704             (bGStatEveryStep ||
1705              (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1706         {
1707             /* Store the pressure in t_state for pressure coupling
1708              * at the next MD step.
1709              */
1710             copy_mat(pres, state->pres_prev);
1711         }
1712
1713         /* #######  END SET VARIABLES FOR NEXT ITERATION ###### */
1714
1715         if ( (membed != NULL) && (!bLastStep) )
1716         {
1717             rescale_membed(step_rel, membed, state_global->x);
1718         }
1719
1720         if (bRerunMD)
1721         {
1722             if (MASTER(cr))
1723             {
1724                 /* read next frame from input trajectory */
1725                 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1726             }
1727
1728             if (PAR(cr))
1729             {
1730                 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1731             }
1732         }
1733
1734         cycles = wallcycle_stop(wcycle, ewcSTEP);
1735         if (DOMAINDECOMP(cr) && wcycle)
1736         {
1737             dd_cycles_add(cr->dd, cycles, ddCyclStep);
1738         }
1739
1740         if (!bRerunMD || !rerun_fr.bStep)
1741         {
1742             /* increase the MD step number */
1743             step++;
1744             step_rel++;
1745         }
1746
1747         /* TODO make a counter-reset module */
1748         /* If it is time to reset counters, set a flag that remains
1749            true until counters actually get reset */
1750         if (step_rel == wcycle_get_reset_counters(wcycle) ||
1751             signals[eglsRESETCOUNTERS].set != 0)
1752         {
1753             if (pme_loadbal_is_active(pme_loadbal))
1754             {
1755                 /* Do not permit counter reset while PME load
1756                  * balancing is active. The only purpose for resetting
1757                  * counters is to measure reliable performance data,
1758                  * and that can't be done before balancing
1759                  * completes.
1760                  *
1761                  * TODO consider fixing this by delaying the reset
1762                  * until after load balancing completes,
1763                  * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1764                 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1765                           "reset mdrun counters at step %" GMX_PRId64 ". Try "
1766                           "resetting counters later in the run, e.g. with gmx "
1767                           "mdrun -resetstep.", step);
1768             }
1769             reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1770                                use_GPU(fr->nbv) ? fr->nbv : NULL);
1771             wcycle_set_reset_counters(wcycle, -1);
1772             if (!(cr->duty & DUTY_PME))
1773             {
1774                 /* Tell our PME node to reset its counters */
1775                 gmx_pme_send_resetcounters(cr, step);
1776             }
1777             /* Correct max_hours for the elapsed time */
1778             max_hours                -= elapsed_time/(60.0*60.0);
1779             /* If mdrun -maxh -resethway was active, it can only trigger once */
1780             bResetCountersHalfMaxH    = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1781             /* Reset can only happen once, so clear the triggering flag. */
1782             signals[eglsRESETCOUNTERS].set = 0;
1783         }
1784
1785         /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1786         IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1787
1788     }
1789     /* End of main MD loop */
1790
1791     /* Closing TNG files can include compressing data. Therefore it is good to do that
1792      * before stopping the time measurements. */
1793     mdoutf_tng_close(outf);
1794
1795     /* Stop measuring walltime */
1796     walltime_accounting_end(walltime_accounting);
1797
1798     if (bRerunMD && MASTER(cr))
1799     {
1800         close_trj(status);
1801     }
1802
1803     if (!(cr->duty & DUTY_PME))
1804     {
1805         /* Tell the PME only node to finish */
1806         gmx_pme_send_finish(cr);
1807     }
1808
1809     if (MASTER(cr))
1810     {
1811         if (ir->nstcalcenergy > 0 && !bRerunMD)
1812         {
1813             print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1814                        eprAVER, mdebin, fcd, groups, &(ir->opts));
1815         }
1816     }
1817
1818     done_mdoutf(outf);
1819
1820     if (bPMETune)
1821     {
1822         pme_loadbal_done(pme_loadbal, cr, fplog, use_GPU(fr->nbv));
1823     }
1824
1825     done_shellfc(fplog, shellfc, step_rel);
1826
1827     if (repl_ex_nst > 0 && MASTER(cr))
1828     {
1829         print_replica_exchange_statistics(fplog, repl_ex);
1830     }
1831
1832     // Clean up swapcoords
1833     if (ir->eSwapCoords != eswapNO)
1834     {
1835         finish_swapcoords(ir->swap);
1836     }
1837
1838     /* IMD cleanup, if bIMD is TRUE. */
1839     IMD_finalize(ir->bIMD, ir->imd);
1840
1841     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1842
1843     return 0;
1844 }