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