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