1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
40 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
83 #include "compute_io.h"
85 #include "checkpoint.h"
86 #include "mtop_util.h"
87 #include "sighandler.h"
103 double do_md(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
104 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
106 gmx_vsite_t *vsite,gmx_constr_t constr,
107 int stepout,t_inputrec *ir,
108 gmx_mtop_t *top_global,
110 t_state *state_global,
112 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
113 gmx_edsam_t ed,t_forcerec *fr,
114 int repl_ex_nst,int repl_ex_seed,gmx_membed_t *membed,
115 real cpt_period,real max_hours,
116 const char *deviceOptions,
118 gmx_runtime_t *runtime)
121 gmx_large_int_t step,step_rel;
124 gmx_bool bGStatEveryStep,bGStat,bNstEner,bCalcEnerPres;
125 gmx_bool bNS,bNStList,bSimAnn,bStopCM,bRerunMD,bNotLastFrame=FALSE,
126 bFirstStep,bStateFromTPX,bInitStep,bLastStep,
127 bBornRadii,bStartingFromCpt;
128 gmx_bool bDoDHDL=FALSE;
129 gmx_bool do_ene,do_log,do_verbose,bRerunWarnNoV=TRUE,
130 bForceUpdate=FALSE,bCPT;
132 gmx_bool bMasterState;
133 int force_flags,cglo_flags;
134 tensor force_vir,shake_vir,total_vir,tmp_vir,pres;
139 t_state *bufstate=NULL;
140 matrix *scale_tot,pcoupl_mu,M,ebox;
143 gmx_repl_ex_t repl_ex=NULL;
147 t_mdebin *mdebin=NULL;
152 gmx_enerdata_t *enerd;
154 gmx_global_stat_t gstat;
155 gmx_update_t upd=NULL;
160 gmx_groups_t *groups;
161 gmx_ekindata_t *ekind, *ekind_save;
162 gmx_shellfc_t shellfc;
163 int count,nconverged=0;
166 gmx_bool bIonize=FALSE;
167 gmx_bool bTCR=FALSE,bConverged=TRUE,bOK,bSumEkinhOld,bExchanged;
169 gmx_bool bResetCountersHalfMaxH=FALSE;
170 gmx_bool bVV,bIterations,bFirstIterate,bTemp,bPres,bTrotter;
171 real temp0,mu_aver=0,dvdl;
173 atom_id *grpindex=NULL;
175 t_coupl_rec *tcr=NULL;
176 rvec *xcopy=NULL,*vcopy=NULL,*cbuf=NULL;
177 matrix boxcopy={{0}},lastbox;
179 real fom,oldfom,veta_save,pcurr,scalevir,tracevir;
182 real saved_conserved_quantity = 0;
187 char sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
188 int handled_stop_condition=gmx_stop_cond_none; /* compare to get_stop_condition*/
189 gmx_iterate_t iterate;
190 gmx_large_int_t multisim_nsteps=-1; /* number of steps to do before first multisim
191 simulation stops. If equal to zero, don't
192 communicate any more between multisims.*/
194 /* Temporary addition for FAHCORE checkpointing */
198 /* Check for special mdrun options */
199 bRerunMD = (Flags & MD_RERUN);
200 bIonize = (Flags & MD_IONIZE);
201 bFFscan = (Flags & MD_FFSCAN);
202 bAppend = (Flags & MD_APPENDFILES);
203 if (Flags & MD_RESETCOUNTERSHALFWAY)
207 /* Signal to reset the counters half the simulation steps. */
208 wcycle_set_reset_counters(wcycle,ir->nsteps/2);
210 /* Signal to reset the counters halfway the simulation time. */
211 bResetCountersHalfMaxH = (max_hours > 0);
214 /* md-vv uses averaged full step velocities for T-control
215 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
216 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
218 if (bVV) /* to store the initial velocities while computing virial */
220 snew(cbuf,top_global->natoms);
222 /* all the iteratative cases - only if there are constraints */
223 bIterations = ((IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
224 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || (IR_NVT_TROTTER(ir))));
228 /* Since we don't know if the frames read are related in any way,
229 * rebuild the neighborlist at every step.
232 ir->nstcalcenergy = 1;
236 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
238 nstglobalcomm = check_nstglobalcomm(fplog,cr,nstglobalcomm,ir);
239 bGStatEveryStep = (nstglobalcomm == 1);
241 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
244 "To reduce the energy communication with nstlist = -1\n"
245 "the neighbor list validity should not be checked at every step,\n"
246 "this means that exact integration is not guaranteed.\n"
247 "The neighbor list validity is checked after:\n"
248 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
249 "In most cases this will result in exact integration.\n"
250 "This reduces the energy communication by a factor of 2 to 3.\n"
251 "If you want less energy communication, set nstlist > 3.\n\n");
254 if (bRerunMD || bFFscan)
258 groups = &top_global->groups;
261 init_md(fplog,cr,ir,oenv,&t,&t0,&state_global->lambda,&lam0,
262 nrnb,top_global,&upd,
263 nfile,fnm,&outf,&mdebin,
264 force_vir,shake_vir,mu_tot,&bSimAnn,&vcm,state_global,Flags);
266 clear_mat(total_vir);
268 /* Energy terms and groups */
270 init_enerdata(top_global->groups.grps[egcENER].nr,ir->n_flambda,enerd);
271 if (DOMAINDECOMP(cr))
277 snew(f,top_global->natoms);
280 /* Kinetic energy data */
282 init_ekindata(fplog,top_global,&(ir->opts),ekind);
283 /* needed for iteration of constraints */
285 init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
286 /* Copy the cos acceleration to the groups struct */
287 ekind->cosacc.cos_accel = ir->cos_accel;
289 gstat = global_stat_init(ir);
292 /* Check for polarizable models and flexible constraints */
293 shellfc = init_shell_flexcon(fplog,
294 top_global,n_flexible_constraints(constr),
295 (ir->bContinuation ||
296 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
297 NULL : state_global->x);
302 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
304 set_deform_reference_box(upd,
305 deform_init_init_step_tpx,
306 deform_init_box_tpx);
308 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
313 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
314 if ((io > 2000) && MASTER(cr))
316 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
320 if (DOMAINDECOMP(cr)) {
321 top = dd_init_local_top(top_global);
324 dd_init_local_state(cr->dd,state_global,state);
326 if (DDMASTER(cr->dd) && ir->nstfout) {
327 snew(f_global,state_global->natoms);
331 /* Initialize the particle decomposition and split the topology */
332 top = split_system(fplog,top_global,ir,cr);
334 pd_cg_range(cr,&fr->cg0,&fr->hcg);
335 pd_at_range(cr,&a0,&a1);
337 top = gmx_mtop_generate_local_top(top_global,ir);
340 a1 = top_global->natoms;
343 state = partdec_init_local_state(cr,state_global);
346 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
349 set_vsite_top(vsite,top,mdatoms,cr);
352 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols) {
353 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
357 make_local_shells(cr,mdatoms,shellfc);
360 if (ir->pull && PAR(cr)) {
361 dd_make_local_pull_groups(NULL,ir->pull,mdatoms);
365 if (DOMAINDECOMP(cr))
367 /* Distribute the charge groups over the nodes from the master node */
368 dd_partition_system(fplog,ir->init_step,cr,TRUE,1,
369 state_global,top_global,ir,
370 state,&f,mdatoms,top,fr,
371 vsite,shellfc,constr,
375 update_mdatoms(mdatoms,state->lambda);
379 if (opt2bSet("-cpi",nfile,fnm))
381 /* Update mdebin with energy history if appending to output files */
382 if ( Flags & MD_APPENDFILES )
384 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
388 /* We might have read an energy history from checkpoint,
389 * free the allocated memory and reset the counts.
391 done_energyhistory(&state_global->enerhist);
392 init_energyhistory(&state_global->enerhist);
395 /* Set the initial energy history in state by updating once */
396 update_energyhistory(&state_global->enerhist,mdebin);
399 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG)) {
400 /* Set the random state if we read a checkpoint file */
401 set_stochd_state(upd,state);
404 /* Initialize constraints */
406 if (!DOMAINDECOMP(cr))
407 set_constraints(constr,top,ir,mdatoms,cr);
410 /* Check whether we have to GCT stuff */
411 bTCR = ftp2bSet(efGCT,nfile,fnm);
414 fprintf(stderr,"Will do General Coupling Theory!\n");
416 gnx = top_global->mols.nr;
418 for(i=0; (i<gnx); i++) {
425 /* We need to be sure replica exchange can only occur
426 * when the energies are current */
427 check_nst_param(fplog,cr,"nstcalcenergy",ir->nstcalcenergy,
428 "repl_ex_nst",&repl_ex_nst);
429 /* This check needs to happen before inter-simulation
430 * signals are initialized, too */
432 if (repl_ex_nst > 0 && MASTER(cr))
433 repl_ex = init_replica_exchange(fplog,cr->ms,state_global,ir,
434 repl_ex_nst,repl_ex_seed);
436 if (!ir->bContinuation && !bRerunMD)
438 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
440 /* Set the velocities of frozen particles to zero */
441 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
445 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
455 /* Constrain the initial coordinates and velocities */
456 do_constrain_first(fplog,constr,ir,mdatoms,state,f,
457 graph,cr,nrnb,fr,top,shake_vir);
461 /* Construct the virtual sites for the initial configuration */
462 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
463 top->idef.iparams,top->idef.il,
464 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
470 /* I'm assuming we need global communication the first time! MRS */
471 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
472 | (bVV ? CGLO_PRESSURE:0)
473 | (bVV ? CGLO_CONSTRAINT:0)
474 | (bRerunMD ? CGLO_RERUNMD:0)
475 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN:0));
477 bSumEkinhOld = FALSE;
478 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
479 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
480 constr,NULL,FALSE,state->box,
481 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,cglo_flags);
482 if (ir->eI == eiVVAK) {
483 /* a second call to get the half step temperature initialized as well */
484 /* we do the same call as above, but turn the pressure off -- internally to
485 compute_globals, this is recognized as a velocity verlet half-step
486 kinetic energy calculation. This minimized excess variables, but
487 perhaps loses some logic?*/
489 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
490 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
491 constr,NULL,FALSE,state->box,
492 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
493 cglo_flags &~ CGLO_PRESSURE);
496 /* Calculate the initial half step temperature, and save the ekinh_old */
497 if (!(Flags & MD_STARTFROMCPT))
499 for(i=0; (i<ir->opts.ngtc); i++)
501 copy_mat(ekind->tcstat[i].ekinh,ekind->tcstat[i].ekinh_old);
506 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
507 and there is no previous step */
509 temp0 = enerd->term[F_TEMP];
511 /* if using an iterative algorithm, we need to create a working directory for the state. */
514 bufstate = init_bufstate(state);
518 snew(xcopy,state->natoms);
519 snew(vcopy,state->natoms);
520 copy_rvecn(state->x,xcopy,0,state->natoms);
521 copy_rvecn(state->v,vcopy,0,state->natoms);
522 copy_mat(state->box,boxcopy);
525 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
526 temperature control */
527 trotter_seq = init_npt_vars(ir,state,&MassQ,bTrotter);
531 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
534 "RMS relative constraint deviation after constraining: %.2e\n",
535 constr_rmsd(constr,FALSE));
537 fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
540 fprintf(stderr,"starting md rerun '%s', reading coordinates from"
541 " input trajectory '%s'\n\n",
542 *(top_global->name),opt2fn("-rerun",nfile,fnm));
545 fprintf(stderr,"Calculated time to finish depends on nsteps from "
546 "run input file,\nwhich may not correspond to the time "
547 "needed to process input trajectory.\n\n");
553 fprintf(stderr,"starting mdrun '%s'\n",
554 *(top_global->name));
557 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
561 sprintf(tbuf,"%s","infinite");
563 if (ir->init_step > 0)
565 fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
566 gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
567 gmx_step_str(ir->init_step,sbuf2),
568 ir->init_step*ir->delta_t);
572 fprintf(stderr,"%s steps, %s ps.\n",
573 gmx_step_str(ir->nsteps,sbuf),tbuf);
579 /* Set and write start time */
580 runtime_start(runtime);
581 print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
582 wallcycle_start(wcycle,ewcRUN);
586 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
588 chkpt_ret=fcCheckPointParallel( cr->nodeid,
590 if ( chkpt_ret == 0 )
591 gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", 0 );
595 /***********************************************************
599 ************************************************************/
601 /* if rerunMD then read coordinates and velocities from input trajectory */
604 if (getenv("GMX_FORCE_UPDATE"))
612 bNotLastFrame = read_first_frame(oenv,&status,
613 opt2fn("-rerun",nfile,fnm),
614 &rerun_fr,TRX_NEED_X | TRX_READ_V);
615 if (rerun_fr.natoms != top_global->natoms)
618 "Number of atoms in trajectory (%d) does not match the "
619 "run input file (%d)\n",
620 rerun_fr.natoms,top_global->natoms);
622 if (ir->ePBC != epbcNONE)
626 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);
628 if (max_cutoff2(ir->ePBC,rerun_fr.box) < sqr(fr->rlistlong))
630 gmx_fatal(FARGS,"Rerun trajectory frame step %d time %f has too small box dimensions",rerun_fr.step,rerun_fr.time);
637 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
640 if (ir->ePBC != epbcNONE)
642 /* Set the shift vectors.
643 * Necessary here when have a static box different from the tpr box.
645 calc_shifts(rerun_fr.box,fr->shift_vec);
649 /* loop over MD steps or if rerunMD to end of input trajectory */
651 /* Skip the first Nose-Hoover integration when we get the state from tpx */
652 bStateFromTPX = !opt2bSet("-cpi",nfile,fnm);
653 bInitStep = bFirstStep && (bStateFromTPX || bVV);
654 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
656 bSumEkinhOld = FALSE;
659 init_global_signals(&gs,cr,ir,repl_ex_nst);
661 step = ir->init_step;
664 if (ir->nstlist == -1)
666 init_nlistheuristics(&nlh,bGStatEveryStep,step);
669 if (MULTISIM(cr) && (repl_ex_nst <=0 ))
671 /* check how many steps are left in other sims */
672 multisim_nsteps=get_multisim_nsteps(cr, ir->nsteps);
676 /* and stop now if we should */
677 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
678 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
679 while (!bLastStep || (bRerunMD && bNotLastFrame)) {
681 wallcycle_start(wcycle,ewcSTEP);
684 if (rerun_fr.bStep) {
685 step = rerun_fr.step;
686 step_rel = step - ir->init_step;
688 if (rerun_fr.bTime) {
698 bLastStep = (step_rel == ir->nsteps);
699 t = t0 + step*ir->delta_t;
702 if (ir->efep != efepNO)
704 if (bRerunMD && rerun_fr.bLambda && (ir->delta_lambda!=0))
706 state_global->lambda = rerun_fr.lambda;
710 state_global->lambda = lam0 + step*ir->delta_lambda;
712 state->lambda = state_global->lambda;
713 bDoDHDL = do_per_step(step,ir->nstdhdl);
718 update_annealing_target_temp(&(ir->opts),t);
723 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
725 for(i=0; i<state_global->natoms; i++)
727 copy_rvec(rerun_fr.x[i],state_global->x[i]);
731 for(i=0; i<state_global->natoms; i++)
733 copy_rvec(rerun_fr.v[i],state_global->v[i]);
738 for(i=0; i<state_global->natoms; i++)
740 clear_rvec(state_global->v[i]);
744 fprintf(stderr,"\nWARNING: Some frames do not contain velocities.\n"
745 " Ekin, temperature and pressure are incorrect,\n"
746 " the virial will be incorrect when constraints are present.\n"
748 bRerunWarnNoV = FALSE;
752 copy_mat(rerun_fr.box,state_global->box);
753 copy_mat(state_global->box,state->box);
755 if (vsite && (Flags & MD_RERUN_VSITE))
757 if (DOMAINDECOMP(cr))
759 gmx_fatal(FARGS,"Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
763 /* Following is necessary because the graph may get out of sync
764 * with the coordinates if we only have every N'th coordinate set
766 mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
767 shift_self(graph,state->box,state->x);
769 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
770 top->idef.iparams,top->idef.il,
771 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
774 unshift_self(graph,state->box,state->x);
779 /* Stop Center of Mass motion */
780 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step,ir->nstcomm));
782 /* Copy back starting coordinates in case we're doing a forcefield scan */
785 for(ii=0; (ii<state->natoms); ii++)
787 copy_rvec(xcopy[ii],state->x[ii]);
788 copy_rvec(vcopy[ii],state->v[ii]);
790 copy_mat(boxcopy,state->box);
795 /* for rerun MD always do Neighbour Searching */
796 bNS = (bFirstStep || ir->nstlist != 0);
801 /* Determine whether or not to do Neighbour Searching and LR */
802 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
804 bNS = (bFirstStep || bExchanged || bNStList ||
805 (ir->nstlist == -1 && nlh.nabnsb > 0));
807 if (bNS && ir->nstlist == -1)
809 set_nlistheuristics(&nlh,bFirstStep || bExchanged,step);
813 /* check whether we should stop because another simulation has
817 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
818 (multisim_nsteps != ir->nsteps) )
825 "Stopping simulation %d because another one has finished\n",
829 gs.sig[eglsCHKPT] = 1;
834 /* < 0 means stop at next step, > 0 means stop at next NS step */
835 if ( (gs.set[eglsSTOPCOND] < 0 ) ||
836 ( (gs.set[eglsSTOPCOND] > 0 ) && ( bNS || ir->nstlist==0)) )
841 /* Determine whether or not to update the Born radii if doing GB */
842 bBornRadii=bFirstStep;
843 if (ir->implicit_solvent && (step % ir->nstgbradii==0))
848 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
849 do_verbose = bVerbose &&
850 (step % stepout == 0 || bFirstStep || bLastStep);
852 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
860 bMasterState = FALSE;
861 /* Correct the new box if it is too skewed */
862 if (DYNAMIC_BOX(*ir))
864 if (correct_box(fplog,step,state->box,graph))
869 if (DOMAINDECOMP(cr) && bMasterState)
871 dd_collect_state(cr->dd,state,state_global);
875 if (DOMAINDECOMP(cr))
877 /* Repartition the domain decomposition */
878 wallcycle_start(wcycle,ewcDOMDEC);
879 dd_partition_system(fplog,step,cr,
880 bMasterState,nstglobalcomm,
881 state_global,top_global,ir,
882 state,&f,mdatoms,top,fr,
883 vsite,shellfc,constr,
884 nrnb,wcycle,do_verbose);
885 wallcycle_stop(wcycle,ewcDOMDEC);
886 /* If using an iterative integrator, reallocate space to match the decomposition */
890 if (MASTER(cr) && do_log && !bFFscan)
892 print_ebin_header(fplog,step,t,state->lambda);
895 if (ir->efep != efepNO)
897 update_mdatoms(mdatoms,state->lambda);
900 if (bRerunMD && rerun_fr.bV)
903 /* We need the kinetic energy at minus the half step for determining
904 * the full step kinetic energy and possibly for T-coupling.*/
905 /* This may not be quite working correctly yet . . . . */
906 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
907 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
908 constr,NULL,FALSE,state->box,
909 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
910 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
912 clear_mat(force_vir);
914 /* Ionize the atoms if necessary */
917 ionize(fplog,oenv,mdatoms,top_global,t,ir,state->x,state->v,
918 mdatoms->start,mdatoms->start+mdatoms->homenr,state->box,cr);
921 /* Update force field in ffscan program */
924 if (update_forcefield(fplog,
926 mdatoms->nr,state->x,state->box)) {
927 if (gmx_parallel_env_initialized())
935 /* We write a checkpoint at this MD step when:
936 * either at an NS step when we signalled through gs,
937 * or at the last step (but not when we do not want confout),
938 * but never at the first step or with rerun.
940 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
941 (bLastStep && (Flags & MD_CONFOUT))) &&
942 step > ir->init_step && !bRerunMD);
945 gs.set[eglsCHKPT] = 0;
948 /* Determine the energy and pressure:
949 * at nstcalcenergy steps and at energy output steps (set below).
951 bNstEner = do_per_step(step,ir->nstcalcenergy);
954 (ir->epc != epcNO && do_per_step(step,ir->nstpcouple)));
956 /* Do we need global communication ? */
957 bGStat = (bCalcEnerPres || bStopCM ||
958 do_per_step(step,nstglobalcomm) ||
959 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
961 do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
963 if (do_ene || do_log)
965 bCalcEnerPres = TRUE;
969 /* these CGLO_ options remain the same throughout the iteration */
970 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
971 (bStopCM ? CGLO_STOPCM : 0) |
972 (bGStat ? CGLO_GSTAT : 0)
975 force_flags = (GMX_FORCE_STATECHANGED |
976 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
977 GMX_FORCE_ALLFORCES |
978 (bNStList ? GMX_FORCE_DOLR : 0) |
980 (bCalcEnerPres ? GMX_FORCE_VIRIAL : 0) |
981 (bDoDHDL ? GMX_FORCE_DHDL : 0)
986 /* Now is the time to relax the shells */
987 count=relax_shell_flexcon(fplog,cr,bVerbose,bFFscan ? step+1 : step,
989 bStopCM,top,top_global,
991 state,f,force_vir,mdatoms,
992 nrnb,wcycle,graph,groups,
993 shellfc,fr,bBornRadii,t,mu_tot,
994 state->natoms,&bConverged,vsite,
1005 /* The coordinates (x) are shifted (to get whole molecules)
1007 * This is parallellized as well, and does communication too.
1008 * Check comments in sim_util.c
1011 do_force(fplog,cr,ir,step,nrnb,wcycle,top,top_global,groups,
1012 state->box,state->x,&state->hist,
1013 f,force_vir,mdatoms,enerd,fcd,
1014 state->lambda,graph,
1015 fr,vsite,mu_tot,t,outf->fp_field,ed,bBornRadii,
1016 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1021 mu_aver = calc_mu_aver(cr,state->x,mdatoms->chargeA,
1022 mu_tot,&top_global->mols,mdatoms,gnx,grpindex);
1025 if (bTCR && bFirstStep)
1027 tcr=init_coupling(fplog,nfile,fnm,cr,fr,mdatoms,&(top->idef));
1028 fprintf(fplog,"Done init_coupling\n");
1032 if (bVV && !bStartingFromCpt && !bRerunMD)
1033 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1035 if (ir->eI==eiVV && bInitStep)
1037 /* if using velocity verlet with full time step Ekin,
1038 * take the first half step only to compute the
1039 * virial for the first step. From there,
1040 * revert back to the initial coordinates
1041 * so that the input is actually the initial step.
1043 copy_rvecn(state->v,cbuf,0,state->natoms); /* should make this better for parallelizing? */
1045 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1046 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ1);
1049 update_coords(fplog,step,ir,mdatoms,state,
1050 f,fr->bTwinRange && bNStList,fr->f_twin,fcd,
1051 ekind,M,wcycle,upd,bInitStep,etrtVELOCITY1,
1052 cr,nrnb,constr,&top->idef);
1056 gmx_iterate_init(&iterate,bIterations && !bInitStep);
1058 /* for iterations, we save these vectors, as we will be self-consistently iterating
1061 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1063 /* save the state */
1064 if (bIterations && iterate.bIterate) {
1065 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1068 bFirstIterate = TRUE;
1069 while (bFirstIterate || (bIterations && iterate.bIterate))
1071 if (bIterations && iterate.bIterate)
1073 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1074 if (bFirstIterate && bTrotter)
1076 /* The first time through, we need a decent first estimate
1077 of veta(t+dt) to compute the constraints. Do
1078 this by computing the box volume part of the
1079 trotter integration at this time. Nothing else
1080 should be changed by this routine here. If
1081 !(first time), we start with the previous value
1084 veta_save = state->veta;
1085 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ0);
1086 vetanew = state->veta;
1087 state->veta = veta_save;
1092 if ( !bRerunMD || rerun_fr.bV || bForceUpdate) { /* Why is rerun_fr.bV here? Unclear. */
1095 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,state,graph,f,
1096 &top->idef,shake_vir,NULL,
1097 cr,nrnb,wcycle,upd,constr,
1098 bInitStep,TRUE,bCalcEnerPres,vetanew);
1100 if (!bOK && !bFFscan)
1102 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1107 { /* Need to unshift here if a do_force has been
1108 called in the previous step */
1109 unshift_self(graph,state->box,state->x);
1113 /* if VV, compute the pressure and constraints */
1114 /* For VV2, we strictly only need this if using pressure
1115 * control, but we really would like to have accurate pressures
1117 * Think about ways around this in the future?
1118 * For now, keep this choice in comments.
1120 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1121 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1123 bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK));
1124 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1125 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1126 constr,NULL,FALSE,state->box,
1127 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1130 | (bTemp ? CGLO_TEMPERATURE:0)
1131 | (bPres ? CGLO_PRESSURE : 0)
1132 | (bPres ? CGLO_CONSTRAINT : 0)
1133 | ((bIterations && iterate.bIterate) ? CGLO_ITERATE : 0)
1134 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1137 /* explanation of above:
1138 a) We compute Ekin at the full time step
1139 if 1) we are using the AveVel Ekin, and it's not the
1140 initial step, or 2) if we are using AveEkin, but need the full
1141 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1142 b) If we are using EkinAveEkin for the kinetic energy for the temperture control, we still feed in
1143 EkinAveVel because it's needed for the pressure */
1145 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1150 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ2);
1154 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1159 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1160 state->veta,&vetanew))
1164 bFirstIterate = FALSE;
1167 if (bTrotter && !bInitStep) {
1168 copy_mat(shake_vir,state->svir_prev);
1169 copy_mat(force_vir,state->fvir_prev);
1170 if (IR_NVT_TROTTER(ir) && ir->eI==eiVV) {
1171 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1172 enerd->term[F_TEMP] = sum_ekin(&(ir->opts),ekind,NULL,(ir->eI==eiVV),FALSE,FALSE);
1173 enerd->term[F_EKIN] = trace(ekind->ekin);
1176 /* if it's the initial step, we performed this first step just to get the constraint virial */
1177 if (bInitStep && ir->eI==eiVV) {
1178 copy_rvecn(cbuf,state->v,0,state->natoms);
1181 if (fr->bSepDVDL && fplog && do_log)
1183 fprintf(fplog,sepdvdlformat,"Constraint",0.0,dvdl);
1185 enerd->term[F_DHDL_CON] += dvdl;
1188 /* MRS -- now done iterating -- compute the conserved quantity */
1190 saved_conserved_quantity = compute_conserved_from_auxiliary(ir,state,&MassQ);
1193 last_ekin = enerd->term[F_EKIN]; /* does this get preserved through checkpointing? */
1195 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1197 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1201 /* ######## END FIRST UPDATE STEP ############## */
1202 /* ######## If doing VV, we now have v(dt) ###### */
1204 /* ################## START TRAJECTORY OUTPUT ################# */
1206 /* Now we have the energies and forces corresponding to the
1207 * coordinates at time t. We must output all of this before
1209 * for RerunMD t is read from input trajectory
1212 if (do_per_step(step,ir->nstxout)) { mdof_flags |= MDOF_X; }
1213 if (do_per_step(step,ir->nstvout)) { mdof_flags |= MDOF_V; }
1214 if (do_per_step(step,ir->nstfout)) { mdof_flags |= MDOF_F; }
1215 if (do_per_step(step,ir->nstxtcout)) { mdof_flags |= MDOF_XTC; }
1216 if (bCPT) { mdof_flags |= MDOF_CPT; };
1218 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1221 /* Enforce writing positions and velocities at end of run */
1222 mdof_flags |= (MDOF_X | MDOF_V);
1227 fcReportProgress( ir->nsteps, step );
1229 /* sync bCPT and fc record-keeping */
1230 if (bCPT && MASTER(cr))
1231 fcRequestCheckPoint();
1234 if (mdof_flags != 0)
1236 wallcycle_start(wcycle,ewcTRAJ);
1239 if (state->flags & (1<<estLD_RNG))
1241 get_stochd_state(upd,state);
1247 state_global->ekinstate.bUpToDate = FALSE;
1251 update_ekinstate(&state_global->ekinstate,ekind);
1252 state_global->ekinstate.bUpToDate = TRUE;
1254 update_energyhistory(&state_global->enerhist,mdebin);
1257 write_traj(fplog,cr,outf,mdof_flags,top_global,
1258 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
1265 if (bLastStep && step_rel == ir->nsteps &&
1266 (Flags & MD_CONFOUT) && MASTER(cr) &&
1267 !bRerunMD && !bFFscan)
1269 /* x and v have been collected in write_traj,
1270 * because a checkpoint file will always be written
1273 fprintf(stderr,"\nWriting final coordinates.\n");
1274 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols &&
1277 /* Make molecules whole only for confout writing */
1278 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
1280 write_sto_conf_mtop(ftp2fn(efSTO,nfile,fnm),
1281 *top_global->name,top_global,
1282 state_global->x,state_global->v,
1283 ir->ePBC,state->box);
1286 wallcycle_stop(wcycle,ewcTRAJ);
1289 /* kludge -- virial is lost with restart for NPT control. Must restart */
1290 if (bStartingFromCpt && bVV)
1292 copy_mat(state->svir_prev,shake_vir);
1293 copy_mat(state->fvir_prev,force_vir);
1295 /* ################## END TRAJECTORY OUTPUT ################ */
1297 /* Determine the wallclock run time up till now */
1298 run_time = gmx_gettime() - (double)runtime->real;
1300 /* Check whether everything is still allright */
1301 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1307 /* this is just make gs.sig compatible with the hack
1308 of sending signals around by MPI_Reduce with together with
1310 if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
1311 gs.sig[eglsSTOPCOND]=1;
1312 if ( gmx_get_stop_condition() == gmx_stop_cond_next )
1313 gs.sig[eglsSTOPCOND]=-1;
1314 /* < 0 means stop at next step, > 0 means stop at next NS step */
1318 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1319 gmx_get_signal_name(),
1320 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1324 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1325 gmx_get_signal_name(),
1326 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1328 handled_stop_condition=(int)gmx_get_stop_condition();
1330 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1331 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1332 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1334 /* Signal to terminate the run */
1335 gs.sig[eglsSTOPCOND] = 1;
1338 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1340 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1343 if (bResetCountersHalfMaxH && MASTER(cr) &&
1344 run_time > max_hours*60.0*60.0*0.495)
1346 gs.sig[eglsRESETCOUNTERS] = 1;
1349 if (ir->nstlist == -1 && !bRerunMD)
1351 /* When bGStatEveryStep=FALSE, global_stat is only called
1352 * when we check the atom displacements, not at NS steps.
1353 * This means that also the bonded interaction count check is not
1354 * performed immediately after NS. Therefore a few MD steps could
1355 * be performed with missing interactions.
1356 * But wrong energies are never written to file,
1357 * since energies are only written after global_stat
1360 if (step >= nlh.step_nscheck)
1362 nlh.nabnsb = natoms_beyond_ns_buffer(ir,fr,&top->cgs,
1363 nlh.scale_tot,state->x);
1367 /* This is not necessarily true,
1368 * but step_nscheck is determined quite conservatively.
1374 /* In parallel we only have to check for checkpointing in steps
1375 * where we do global communication,
1376 * otherwise the other nodes don't know.
1378 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1381 run_time >= nchkpt*cpt_period*60.0)) &&
1382 gs.set[eglsCHKPT] == 0)
1384 gs.sig[eglsCHKPT] = 1;
1389 gmx_iterate_init(&iterate,bIterations);
1392 /* for iterations, we save these vectors, as we will be redoing the calculations */
1393 if (bIterations && iterate.bIterate)
1395 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1397 bFirstIterate = TRUE;
1398 while (bFirstIterate || (bIterations && iterate.bIterate))
1400 /* We now restore these vectors to redo the calculation with improved extended variables */
1403 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1406 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1407 so scroll down for that logic */
1409 /* ######### START SECOND UPDATE STEP ################# */
1410 /* Box is changed in update() when we do pressure coupling,
1411 * but we should still use the old box for energy corrections and when
1412 * writing it to the energy file, so it matches the trajectory files for
1413 * the same timestep above. Make a copy in a separate array.
1415 copy_mat(state->box,lastbox);
1418 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1420 wallcycle_start(wcycle,ewcUPDATE);
1422 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1425 if (bIterations && iterate.bIterate)
1433 /* we use a new value of scalevir to converge the iterations faster */
1434 scalevir = tracevir/trace(shake_vir);
1436 msmul(shake_vir,scalevir,shake_vir);
1437 m_add(force_vir,shake_vir,total_vir);
1438 clear_mat(shake_vir);
1440 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ3);
1441 /* We can only do Berendsen coupling after we have summed
1442 * the kinetic energy or virial. Since the happens
1443 * in global_state after update, we should only do it at
1444 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1449 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1450 update_pcouple(fplog,step,ir,state,pcoupl_mu,M,wcycle,
1456 /* velocity half-step update */
1457 update_coords(fplog,step,ir,mdatoms,state,f,
1458 fr->bTwinRange && bNStList,fr->f_twin,fcd,
1459 ekind,M,wcycle,upd,FALSE,etrtVELOCITY2,
1460 cr,nrnb,constr,&top->idef);
1463 /* Above, initialize just copies ekinh into ekin,
1464 * it doesn't copy position (for VV),
1465 * and entire integrator for MD.
1470 copy_rvecn(state->x,cbuf,0,state->natoms);
1473 update_coords(fplog,step,ir,mdatoms,state,f,fr->bTwinRange && bNStList,fr->f_twin,fcd,
1474 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1475 wallcycle_stop(wcycle,ewcUPDATE);
1477 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,state,graph,f,
1478 &top->idef,shake_vir,force_vir,
1479 cr,nrnb,wcycle,upd,constr,
1480 bInitStep,FALSE,bCalcEnerPres,state->veta);
1484 /* erase F_EKIN and F_TEMP here? */
1485 /* just compute the kinetic energy at the half step to perform a trotter step */
1486 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1487 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1488 constr,NULL,FALSE,lastbox,
1489 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1490 cglo_flags | CGLO_TEMPERATURE
1492 wallcycle_start(wcycle,ewcUPDATE);
1493 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ4);
1494 /* now we know the scaling, we can compute the positions again again */
1495 copy_rvecn(cbuf,state->x,0,state->natoms);
1497 update_coords(fplog,step,ir,mdatoms,state,f,fr->bTwinRange && bNStList,fr->f_twin,fcd,
1498 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1499 wallcycle_stop(wcycle,ewcUPDATE);
1501 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1502 /* are the small terms in the shake_vir here due
1503 * to numerical errors, or are they important
1504 * physically? I'm thinking they are just errors, but not completely sure.
1505 * For now, will call without actually constraining, constr=NULL*/
1506 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,state,graph,f,
1507 &top->idef,tmp_vir,force_vir,
1508 cr,nrnb,wcycle,upd,NULL,
1509 bInitStep,FALSE,bCalcEnerPres,
1512 if (!bOK && !bFFscan)
1514 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1517 if (fr->bSepDVDL && fplog && do_log)
1519 fprintf(fplog,sepdvdlformat,"Constraint",0.0,dvdl);
1521 enerd->term[F_DHDL_CON] += dvdl;
1525 /* Need to unshift here */
1526 unshift_self(graph,state->box,state->x);
1531 wallcycle_start(wcycle,ewcVSITECONSTR);
1534 shift_self(graph,state->box,state->x);
1536 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
1537 top->idef.iparams,top->idef.il,
1538 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
1542 unshift_self(graph,state->box,state->x);
1544 wallcycle_stop(wcycle,ewcVSITECONSTR);
1547 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1548 if (ir->nstlist == -1 && bFirstIterate)
1550 gs.sig[eglsNABNSB] = nlh.nabnsb;
1552 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1553 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1555 bFirstIterate ? &gs : NULL,
1556 (step_rel % gs.nstms == 0) &&
1557 (multisim_nsteps<0 || (step_rel<multisim_nsteps)),
1559 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1561 | (!EI_VV(ir->eI) ? CGLO_ENERGY : 0)
1562 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1563 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1564 | (bIterations && iterate.bIterate ? CGLO_ITERATE : 0)
1565 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1568 if (ir->nstlist == -1 && bFirstIterate)
1570 nlh.nabnsb = gs.set[eglsNABNSB];
1571 gs.set[eglsNABNSB] = 0;
1573 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1574 /* ############# END CALC EKIN AND PRESSURE ################# */
1576 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1577 the virial that should probably be addressed eventually. state->veta has better properies,
1578 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1579 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1582 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1583 trace(shake_vir),&tracevir))
1587 bFirstIterate = FALSE;
1590 update_box(fplog,step,ir,mdatoms,state,graph,f,
1591 ir->nstlist==-1 ? &nlh.scale_tot : NULL,pcoupl_mu,nrnb,wcycle,upd,bInitStep,FALSE);
1593 /* ################# END UPDATE STEP 2 ################# */
1594 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1596 /* The coordinates (x) were unshifted in update */
1597 if (bFFscan && (shellfc==NULL || bConverged))
1599 if (print_forcefield(fplog,enerd->term,mdatoms->homenr,
1601 &(top_global->mols),mdatoms->massT,pres))
1603 if (gmx_parallel_env_initialized())
1607 fprintf(stderr,"\n");
1613 /* We will not sum ekinh_old,
1614 * so signal that we still have to do it.
1616 bSumEkinhOld = TRUE;
1621 /* Only do GCT when the relaxation of shells (minimization) has converged,
1622 * otherwise we might be coupling to bogus energies.
1623 * In parallel we must always do this, because the other sims might
1627 /* Since this is called with the new coordinates state->x, I assume
1628 * we want the new box state->box too. / EL 20040121
1630 do_coupling(fplog,oenv,nfile,fnm,tcr,t,step,enerd->term,fr,
1632 mdatoms,&(top->idef),mu_aver,
1633 top_global->mols.nr,cr,
1634 state->box,total_vir,pres,
1635 mu_tot,state->x,f,bConverged);
1639 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1641 /* sum up the foreign energy and dhdl terms */
1642 sum_dhdl(enerd,state->lambda,ir);
1644 /* use the directly determined last velocity, not actually the averaged half steps */
1645 if (bTrotter && ir->eI==eiVV)
1647 enerd->term[F_EKIN] = last_ekin;
1649 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1653 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1657 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir,state,&MassQ);
1659 /* Check for excessively large energies */
1663 real etot_max = 1e200;
1665 real etot_max = 1e30;
1667 if (fabs(enerd->term[F_ETOT]) > etot_max)
1669 fprintf(stderr,"Energy too large (%g), giving up\n",
1670 enerd->term[F_ETOT]);
1673 /* ######### END PREPARING EDR OUTPUT ########### */
1675 /* Time for performance */
1676 if (((step % stepout) == 0) || bLastStep)
1678 runtime_upd_proc(runtime);
1684 gmx_bool do_dr,do_or;
1686 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1690 upd_mdebin(mdebin,bDoDHDL, TRUE,
1691 t,mdatoms->tmass,enerd,state,lastbox,
1692 shake_vir,force_vir,total_vir,pres,
1693 ekind,mu_tot,constr);
1697 upd_mdebin_step(mdebin);
1700 do_dr = do_per_step(step,ir->nstdisreout);
1701 do_or = do_per_step(step,ir->nstorireout);
1703 print_ebin(outf->fp_ene,do_ene,do_dr,do_or,do_log?fplog:NULL,
1705 eprNORMAL,bCompact,mdebin,fcd,groups,&(ir->opts));
1707 if (ir->ePull != epullNO)
1709 pull_print_output(ir->pull,step,t);
1712 if (do_per_step(step,ir->nstlog))
1714 if(fflush(fplog) != 0)
1716 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of quota?");
1722 /* Remaining runtime */
1723 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal() ))
1727 fprintf(stderr,"\n");
1729 print_time(stderr,runtime,step,ir,cr);
1732 /* Replica exchange */
1734 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1735 do_per_step(step,repl_ex_nst))
1737 bExchanged = replica_exchange(fplog,cr,repl_ex,
1738 state_global,enerd->term,
1741 if (bExchanged && DOMAINDECOMP(cr))
1743 dd_partition_system(fplog,step,cr,TRUE,1,
1744 state_global,top_global,ir,
1745 state,&f,mdatoms,top,fr,
1746 vsite,shellfc,constr,
1753 bStartingFromCpt = FALSE;
1755 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1756 /* With all integrators, except VV, we need to retain the pressure
1757 * at the current step for coupling at the next step.
1759 if ((state->flags & (1<<estPRES_PREV)) &&
1761 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1763 /* Store the pressure in t_state for pressure coupling
1764 * at the next MD step.
1766 copy_mat(pres,state->pres_prev);
1769 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1771 if ( (membed!=NULL) && (!bLastStep) )
1772 rescale_membed(step_rel,membed,state_global->x);
1778 /* read next frame from input trajectory */
1779 bNotLastFrame = read_next_frame(oenv,status,&rerun_fr);
1784 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
1788 if (!bRerunMD || !rerun_fr.bStep)
1790 /* increase the MD step number */
1795 cycles = wallcycle_stop(wcycle,ewcSTEP);
1796 if (DOMAINDECOMP(cr) && wcycle)
1798 dd_cycles_add(cr->dd,cycles,ddCyclStep);
1801 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1802 gs.set[eglsRESETCOUNTERS] != 0)
1804 /* Reset all the counters related to performance over the run */
1805 reset_all_counters(fplog,cr,step,&step_rel,ir,wcycle,nrnb,runtime);
1806 wcycle_set_reset_counters(wcycle,-1);
1807 /* Correct max_hours for the elapsed time */
1808 max_hours -= run_time/(60.0*60.0);
1809 bResetCountersHalfMaxH = FALSE;
1810 gs.set[eglsRESETCOUNTERS] = 0;
1814 /* End of main MD loop */
1818 runtime_end(runtime);
1820 if (bRerunMD && MASTER(cr))
1825 if (!(cr->duty & DUTY_PME))
1827 /* Tell the PME only node to finish */
1833 if (ir->nstcalcenergy > 0 && !bRerunMD)
1835 print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
1836 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
1844 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1846 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)));
1847 fprintf(fplog,"Average number of atoms that crossed the half buffer length: %.1f\n\n",nlh.ab/nlh.nns);
1850 if (shellfc && fplog)
1852 fprintf(fplog,"Fraction of iterations that converged: %.2f %%\n",
1853 (nconverged*100.0)/step_rel);
1854 fprintf(fplog,"Average number of force evaluations per MD step: %.2f\n\n",
1858 if (repl_ex_nst > 0 && MASTER(cr))
1860 print_replica_exchange_statistics(fplog,repl_ex);
1863 runtime->nsteps_done = step_rel;