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