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
70 #include "mpelogging.h"
77 #include "compute_io.h"
79 #include "checkpoint.h"
80 #include "mtop_util.h"
81 #include "sighandler.h"
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,
100 gmx_vsite_t *vsite,gmx_constr_t constr,
101 int stepout,t_inputrec *ir,
102 gmx_mtop_t *top_global,
104 t_state *state_global,
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,
112 gmx_runtime_t *runtime)
115 gmx_large_int_t step,step_rel;
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;
126 gmx_bool bMasterState;
127 int force_flags,cglo_flags;
128 tensor force_vir,shake_vir,total_vir,tmp_vir,pres;
133 t_state *bufstate=NULL;
134 matrix *scale_tot,pcoupl_mu,M,ebox;
137 gmx_repl_ex_t repl_ex=NULL;
140 t_mdebin *mdebin=NULL;
141 df_history_t df_history;
146 gmx_enerdata_t *enerd;
148 gmx_global_stat_t gstat;
149 gmx_update_t upd=NULL;
152 gmx_rng_t mcrng=NULL;
154 gmx_groups_t *groups;
155 gmx_ekindata_t *ekind, *ekind_save;
156 gmx_shellfc_t shellfc;
157 int count,nconverged=0;
160 gmx_bool bIonize=FALSE;
161 gmx_bool bTCR=FALSE,bConverged=TRUE,bOK,bSumEkinhOld,bExchanged;
163 gmx_bool bResetCountersHalfMaxH=FALSE;
164 gmx_bool bVV,bIterations,bFirstIterate,bTemp,bPres,bTrotter;
167 atom_id *grpindex=NULL;
169 t_coupl_rec *tcr=NULL;
170 rvec *xcopy=NULL,*vcopy=NULL,*cbuf=NULL;
171 matrix boxcopy={{0}},lastbox;
173 real fom,oldfom,veta_save,pcurr,scalevir,tracevir;
181 real saved_conserved_quantity = 0;
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.*/
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");
204 /* Temporary addition for FAHCORE checkpointing */
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)
217 /* Signal to reset the counters half the simulation steps. */
218 wcycle_set_reset_counters(wcycle,ir->nsteps/2);
220 /* Signal to reset the counters halfway the simulation time. */
221 bResetCountersHalfMaxH = (max_hours > 0);
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; */
228 if (bVV) /* to store the initial velocities while computing virial */
230 snew(cbuf,top_global->natoms);
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)));
238 /* Since we don't know if the frames read are related in any way,
239 * rebuild the neighborlist at every step.
242 ir->nstcalcenergy = 1;
246 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
248 nstglobalcomm = check_nstglobalcomm(fplog,cr,nstglobalcomm,ir);
249 bGStatEveryStep = (nstglobalcomm == 1);
251 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
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");
264 if (bRerunMD || bFFscan)
268 groups = &top_global->groups;
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);
277 clear_mat(total_vir);
279 /* Energy terms and groups */
281 init_enerdata(top_global->groups.grps[egcENER].nr,ir->fepvals->n_lambda,
283 if (DOMAINDECOMP(cr))
289 snew(f,top_global->natoms);
292 /* lambda Monte carlo random number generator */
295 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
297 /* copy the state into df_history */
298 copy_df_history(&df_history,&state_global->dfhist);
300 /* Kinetic energy data */
302 init_ekindata(fplog,top_global,&(ir->opts),ekind);
303 /* needed for iteration of constraints */
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;
309 gstat = global_stat_init(ir);
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);
321 #ifdef GMX_THREAD_MPI
322 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
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);
333 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
334 if ((io > 2000) && MASTER(cr))
336 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
340 if (DOMAINDECOMP(cr)) {
341 top = dd_init_local_top(top_global);
344 dd_init_local_state(cr->dd,state_global,state);
346 if (DDMASTER(cr->dd) && ir->nstfout) {
347 snew(f_global,state_global->natoms);
351 /* Initialize the particle decomposition and split the topology */
352 top = split_system(fplog,top_global,ir,cr);
354 pd_cg_range(cr,&fr->cg0,&fr->hcg);
355 pd_at_range(cr,&a0,&a1);
357 top = gmx_mtop_generate_local_top(top_global,ir);
360 a1 = top_global->natoms;
363 state = partdec_init_local_state(cr,state_global);
366 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
369 set_vsite_top(vsite,top,mdatoms,cr);
372 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols) {
373 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
377 make_local_shells(cr,mdatoms,shellfc);
380 if (ir->pull && PAR(cr)) {
381 dd_make_local_pull_groups(NULL,ir->pull,mdatoms);
385 if (DOMAINDECOMP(cr))
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,
395 update_mdatoms(mdatoms,state->lambda[efptMASS]);
397 if (opt2bSet("-cpi",nfile,fnm))
399 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr);
403 bStateFromCP = FALSE;
410 /* Update mdebin with energy history if appending to output files */
411 if ( Flags & MD_APPENDFILES )
413 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
417 /* We might have read an energy history from checkpoint,
418 * free the allocated memory and reset the counts.
420 done_energyhistory(&state_global->enerhist);
421 init_energyhistory(&state_global->enerhist);
424 /* Set the initial energy history in state by updating once */
425 update_energyhistory(&state_global->enerhist,mdebin);
428 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
430 /* Set the random state if we read a checkpoint file */
431 set_stochd_state(upd,state);
434 if (state->flags & (1<<estMC_RNG))
436 set_mc_state(mcrng,state);
439 /* Initialize constraints */
441 if (!DOMAINDECOMP(cr))
442 set_constraints(constr,top,ir,mdatoms,cr);
445 /* Check whether we have to GCT stuff */
446 bTCR = ftp2bSet(efGCT,nfile,fnm);
449 fprintf(stderr,"Will do General Coupling Theory!\n");
451 gnx = top_global->mols.nr;
453 for(i=0; (i<gnx); i++) {
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 */
467 if (repl_ex_nst > 0 && MASTER(cr))
469 repl_ex = init_replica_exchange(fplog,cr->ms,state_global,ir,
470 repl_ex_nst,repl_ex_nex,repl_ex_seed);
472 if (!ir->bContinuation && !bRerunMD)
474 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
476 /* Set the velocities of frozen particles to zero */
477 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
481 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
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);
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);
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))
510 nstfep = ir->expandedvals->nstexpanded;
512 if (repl_ex_nst > 0 && repl_ex_nst > nstfep)
514 nstfep = repl_ex_nst;
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));
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?*/
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));
544 /* Calculate the initial half step temperature, and save the ekinh_old */
545 if (!(Flags & MD_STARTFROMCPT))
547 for(i=0; (i<ir->opts.ngtc); i++)
549 copy_mat(ekind->tcstat[i].ekinh,ekind->tcstat[i].ekinh_old);
554 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
555 and there is no previous step */
558 /* if using an iterative algorithm, we need to create a working directory for the state. */
561 bufstate = init_bufstate(state);
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);
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);
578 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
581 "RMS relative constraint deviation after constraining: %.2e\n",
582 constr_rmsd(constr,FALSE));
584 if (EI_STATE_VELOCITY(ir->eI))
586 fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
590 fprintf(stderr,"starting md rerun '%s', reading coordinates from"
591 " input trajectory '%s'\n\n",
592 *(top_global->name),opt2fn("-rerun",nfile,fnm));
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");
603 fprintf(stderr,"starting mdrun '%s'\n",
604 *(top_global->name));
607 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
611 sprintf(tbuf,"%s","infinite");
613 if (ir->init_step > 0)
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);
622 fprintf(stderr,"%s steps, %s ps.\n",
623 gmx_step_str(ir->nsteps,sbuf),tbuf);
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);
638 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
640 chkpt_ret=fcCheckPointParallel( cr->nodeid,
642 if ( chkpt_ret == 0 )
643 gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", 0 );
647 /***********************************************************
651 ************************************************************/
653 /* if rerunMD then read coordinates and velocities from input trajectory */
656 if (getenv("GMX_FORCE_UPDATE"))
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)
670 "Number of atoms in trajectory (%d) does not match the "
671 "run input file (%d)\n",
672 rerun_fr.natoms,top_global->natoms);
674 if (ir->ePBC != epbcNONE)
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);
680 if (max_cutoff2(ir->ePBC,rerun_fr.box) < sqr(fr->rlistlong))
682 gmx_fatal(FARGS,"Rerun trajectory frame step %d time %f has too small box dimensions",rerun_fr.step,rerun_fr.time);
689 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
692 if (ir->ePBC != epbcNONE)
694 /* Set the shift vectors.
695 * Necessary here when have a static box different from the tpr box.
697 calc_shifts(rerun_fr.box,fr->shift_vec);
701 /* loop over MD steps or if rerunMD to end of input trajectory */
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;
708 bSumEkinhOld = FALSE;
711 init_global_signals(&gs,cr,ir,repl_ex_nst);
713 step = ir->init_step;
716 if (ir->nstlist == -1)
718 init_nlistheuristics(&nlh,bGStatEveryStep,step);
721 if (MULTISIM(cr) && (repl_ex_nst <=0 ))
723 /* check how many steps are left in other sims */
724 multisim_nsteps=get_multisim_nsteps(cr, ir->nsteps);
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)) {
733 wallcycle_start(wcycle,ewcSTEP);
735 GMX_MPE_LOG(ev_timestep1);
738 if (rerun_fr.bStep) {
739 step = rerun_fr.step;
740 step_rel = step - ir->init_step;
742 if (rerun_fr.bTime) {
752 bLastStep = (step_rel == ir->nsteps);
753 t = t0 + step*ir->delta_t;
756 if (ir->efep != efepNO || ir->bSimTemp)
758 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
759 requiring different logic. */
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));
769 update_annealing_target_temp(&(ir->opts),t);
774 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
776 for(i=0; i<state_global->natoms; i++)
778 copy_rvec(rerun_fr.x[i],state_global->x[i]);
782 for(i=0; i<state_global->natoms; i++)
784 copy_rvec(rerun_fr.v[i],state_global->v[i]);
789 for(i=0; i<state_global->natoms; i++)
791 clear_rvec(state_global->v[i]);
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"
799 bRerunWarnNoV = FALSE;
803 copy_mat(rerun_fr.box,state_global->box);
804 copy_mat(state_global->box,state->box);
806 if (vsite && (Flags & MD_RERUN_VSITE))
808 if (DOMAINDECOMP(cr))
810 gmx_fatal(FARGS,"Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
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
817 mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
818 shift_self(graph,state->box,state->x);
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);
825 unshift_self(graph,state->box,state->x);
830 /* Stop Center of Mass motion */
831 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step,ir->nstcomm));
833 /* Copy back starting coordinates in case we're doing a forcefield scan */
836 for(ii=0; (ii<state->natoms); ii++)
838 copy_rvec(xcopy[ii],state->x[ii]);
839 copy_rvec(vcopy[ii],state->v[ii]);
841 copy_mat(boxcopy,state->box);
846 /* for rerun MD always do Neighbour Searching */
847 bNS = (bFirstStep || ir->nstlist != 0);
852 /* Determine whether or not to do Neighbour Searching and LR */
853 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
855 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
856 (ir->nstlist == -1 && nlh.nabnsb > 0));
858 if (bNS && ir->nstlist == -1)
860 set_nlistheuristics(&nlh,bFirstStep || bExchanged || bDoFEP, step);
864 /* check whether we should stop because another simulation has
868 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
869 (multisim_nsteps != ir->nsteps) )
876 "Stopping simulation %d because another one has finished\n",
880 gs.sig[eglsCHKPT] = 1;
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)) )
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))
899 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
900 do_verbose = bVerbose &&
901 (step % stepout == 0 || bFirstStep || bLastStep);
903 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
911 bMasterState = FALSE;
912 /* Correct the new box if it is too skewed */
913 if (DYNAMIC_BOX(*ir))
915 if (correct_box(fplog,step,state->box,graph))
920 if (DOMAINDECOMP(cr) && bMasterState)
922 dd_collect_state(cr->dd,state,state_global);
926 if (DOMAINDECOMP(cr))
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 */
941 if (MASTER(cr) && do_log && !bFFscan)
943 print_ebin_header(fplog,step,t,state->lambda[efptFEP]); /* can we improve the information printed here? */
946 if (ir->efep != efepNO)
948 update_mdatoms(mdatoms,state->lambda[efptMASS]);
951 if (bRerunMD && rerun_fr.bV)
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);
963 clear_mat(force_vir);
965 /* Ionize the atoms if necessary */
968 ionize(fplog,oenv,mdatoms,top_global,t,ir,state->x,state->v,
969 mdatoms->start,mdatoms->start+mdatoms->homenr,state->box,cr);
972 /* Update force field in ffscan program */
975 if (update_forcefield(fplog,
977 mdatoms->nr,state->x,state->box))
985 GMX_MPE_LOG(ev_timestep2);
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.
992 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
993 (bLastStep && (Flags & MD_CONFOUT))) &&
994 step > ir->init_step && !bRerunMD);
997 gs.set[eglsCHKPT] = 0;
1000 /* Determine the energy and pressure:
1001 * at nstcalcenergy steps and at energy output steps (set below).
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);
1007 bNstEner = do_per_step(step,ir->nstcalcenergy);
1011 (ir->epc > epcNO && do_per_step(step,ir->nstpcouple)));
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));
1018 do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
1020 if (do_ene || do_log)
1022 bCalcEnerPres = TRUE;
1026 /* these CGLO_ options remain the same throughout the iteration */
1027 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1028 (bGStat ? CGLO_GSTAT : 0)
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) |
1036 (bCalcEnerPres ? GMX_FORCE_VIRIAL : 0) |
1037 (bDoFEP ? GMX_FORCE_DHDL : 0)
1042 /* Now is the time to relax the shells */
1043 count=relax_shell_flexcon(fplog,cr,bVerbose,bFFscan ? step+1 : step,
1045 bStopCM,top,top_global,
1047 state,f,force_vir,mdatoms,
1048 nrnb,wcycle,graph,groups,
1049 shellfc,fr,bBornRadii,t,mu_tot,
1050 state->natoms,&bConverged,vsite,
1061 /* The coordinates (x) are shifted (to get whole molecules)
1063 * This is parallellized as well, and does communication too.
1064 * Check comments in sim_util.c
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);
1074 GMX_BARRIER(cr->mpi_comm_mygroup);
1078 mu_aver = calc_mu_aver(cr,state->x,mdatoms->chargeA,
1079 mu_tot,&top_global->mols,mdatoms,gnx,grpindex);
1082 if (bTCR && bFirstStep)
1084 tcr=init_coupling(fplog,nfile,fnm,cr,fr,mdatoms,&(top->idef));
1085 fprintf(fplog,"Done init_coupling\n");
1089 if (bVV && !bStartingFromCpt && !bRerunMD)
1090 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1092 if (ir->eI==eiVV && bInitStep)
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.
1100 copy_rvecn(state->v,cbuf,0,state->natoms); /* should make this better for parallelizing? */
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);
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);
1113 gmx_iterate_init(&iterate,bIterations && !bInitStep);
1115 /* for iterations, we save these vectors, as we will be self-consistently iterating
1118 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1120 /* save the state */
1121 if (bIterations && iterate.bIterate) {
1122 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1125 bFirstIterate = TRUE;
1126 while (bFirstIterate || (bIterations && iterate.bIterate))
1128 if (bIterations && iterate.bIterate)
1130 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1131 if (bFirstIterate && bTrotter)
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
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;
1149 if ( !bRerunMD || rerun_fr.bV || bForceUpdate) { /* Why is rerun_fr.bV here? Unclear. */
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);
1157 if (!bOK && !bFFscan)
1159 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
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);
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
1174 * Think about ways around this in the future?
1175 * For now, keep this choice in comments.
1177 /* bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1178 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
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?*/
1183 bSumEkinhOld = TRUE;
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,
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)
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 */
1207 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1212 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ2);
1216 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1221 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1222 state->veta,&vetanew))
1226 bFirstIterate = FALSE;
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);
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);
1244 if (fr->bSepDVDL && fplog && do_log)
1246 fprintf(fplog,sepdvdlformat,"Constraint",0.0,dvdl);
1248 enerd->term[F_DVDL_BONDED] += dvdl;
1250 GMX_MPE_LOG(ev_timestep1);
1253 /* MRS -- now done iterating -- compute the conserved quantity */
1255 saved_conserved_quantity = compute_conserved_from_auxiliary(ir,state,&MassQ);
1258 last_ekin = enerd->term[F_EKIN];
1260 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1262 saved_conserved_quantity -= enerd->term[F_DISPCORR];
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);
1268 /* ######## END FIRST UPDATE STEP ############## */
1269 /* ######## If doing VV, we now have v(dt) ###### */
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. */
1276 lamnew = ExpandedEnsembleDynamics(fplog,ir,enerd,state,&MassQ,&df_history,step,mcrng,state->v,mdatoms);
1278 /* ################## START TRAJECTORY OUTPUT ################# */
1280 /* Now we have the energies and forces corresponding to the
1281 * coordinates at time t. We must output all of this before
1283 * for RerunMD t is read from input trajectory
1285 GMX_MPE_LOG(ev_output_start);
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; };
1294 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1297 /* Enforce writing positions and velocities at end of run */
1298 mdof_flags |= (MDOF_X | MDOF_V);
1303 fcReportProgress( ir->nsteps, step );
1305 /* sync bCPT and fc record-keeping */
1306 if (bCPT && MASTER(cr))
1307 fcRequestCheckPoint();
1310 if (mdof_flags != 0)
1312 wallcycle_start(wcycle,ewcTRAJ);
1315 if (state->flags & (1<<estLD_RNG))
1317 get_stochd_state(upd,state);
1319 if (state->flags & (1<<estMC_RNG))
1321 get_mc_state(mcrng,state);
1327 state_global->ekinstate.bUpToDate = FALSE;
1331 update_ekinstate(&state_global->ekinstate,ekind);
1332 state_global->ekinstate.bUpToDate = TRUE;
1334 update_energyhistory(&state_global->enerhist,mdebin);
1335 if (ir->efep!=efepNO || ir->bSimTemp)
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);
1345 write_traj(fplog,cr,outf,mdof_flags,top_global,
1346 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
1353 if (bLastStep && step_rel == ir->nsteps &&
1354 (Flags & MD_CONFOUT) && MASTER(cr) &&
1355 !bRerunMD && !bFFscan)
1357 /* x and v have been collected in write_traj,
1358 * because a checkpoint file will always be written
1361 fprintf(stderr,"\nWriting final coordinates.\n");
1362 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols &&
1365 /* Make molecules whole only for confout writing */
1366 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
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);
1374 wallcycle_stop(wcycle,ewcTRAJ);
1376 GMX_MPE_LOG(ev_output_finish);
1378 /* kludge -- virial is lost with restart for NPT control. Must restart */
1379 if (bStartingFromCpt && bVV)
1381 copy_mat(state->svir_prev,shake_vir);
1382 copy_mat(state->fvir_prev,force_vir);
1384 /* ################## END TRAJECTORY OUTPUT ################ */
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).
1393 bNstEner = (bGStatEveryStep || do_per_step(step,ir->nstcalcenergy));
1394 bCalcEnerPres = bNstEner;
1396 /* Do we need global communication ? */
1397 bGStat = (bGStatEveryStep || bStopCM || bNS ||
1398 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1400 do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
1402 if (do_ene || do_log)
1404 bCalcEnerPres = TRUE;
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
1417 /* this is just make gs.sig compatible with the hack
1418 of sending signals around by MPI_Reduce with together with
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 */
1428 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1429 gmx_get_signal_name(),
1430 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1434 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1435 gmx_get_signal_name(),
1436 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1438 handled_stop_condition=(int)gmx_get_stop_condition();
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)
1444 /* Signal to terminate the run */
1445 gs.sig[eglsSTOPCOND] = 1;
1448 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1450 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1453 if (bResetCountersHalfMaxH && MASTER(cr) &&
1454 run_time > max_hours*60.0*60.0*0.495)
1456 gs.sig[eglsRESETCOUNTERS] = 1;
1459 if (ir->nstlist == -1 && !bRerunMD)
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
1470 if (step >= nlh.step_nscheck)
1472 nlh.nabnsb = natoms_beyond_ns_buffer(ir,fr,&top->cgs,
1473 nlh.scale_tot,state->x);
1477 /* This is not necessarily true,
1478 * but step_nscheck is determined quite conservatively.
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.
1488 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1491 run_time >= nchkpt*cpt_period*60.0)) &&
1492 gs.set[eglsCHKPT] == 0)
1494 gs.sig[eglsCHKPT] = 1;
1498 /* at the start of step, randomize the velocities */
1499 if (ETC_ANDERSEN(ir->etc) && EI_VV(ir->eI))
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)
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);
1515 gmx_iterate_init(&iterate,bIterations);
1518 /* for iterations, we save these vectors, as we will be redoing the calculations */
1519 if (bIterations && iterate.bIterate)
1521 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1523 bFirstIterate = TRUE;
1524 while (bFirstIterate || (bIterations && iterate.bIterate))
1526 /* We now restore these vectors to redo the calculation with improved extended variables */
1529 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1532 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1533 so scroll down for that logic */
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.
1542 copy_mat(state->box,lastbox);
1545 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1547 wallcycle_start(wcycle,ewcUPDATE);
1549 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1552 if (bIterations && iterate.bIterate)
1560 /* we use a new value of scalevir to converge the iterations faster */
1561 scalevir = tracevir/trace(shake_vir);
1563 msmul(shake_vir,scalevir,shake_vir);
1564 m_add(force_vir,shake_vir,total_vir);
1565 clear_mat(shake_vir);
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.
1576 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1577 update_pcouple(fplog,step,ir,state,pcoupl_mu,M,wcycle,
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);
1590 /* Above, initialize just copies ekinh into ekin,
1591 * it doesn't copy position (for VV),
1592 * and entire integrator for MD.
1597 copy_rvecn(state->x,cbuf,0,state->natoms);
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);
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);
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
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);
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);
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,
1639 if (!bOK && !bFFscan)
1641 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1644 if (fr->bSepDVDL && fplog && do_log)
1646 fprintf(fplog,sepdvdlformat,"Constraint dV/dl",0.0,dvdl);
1648 enerd->term[F_DVDL_BONDED] += dvdl;
1652 /* Need to unshift here */
1653 unshift_self(graph,state->box,state->x);
1656 GMX_BARRIER(cr->mpi_comm_mygroup);
1657 GMX_MPE_LOG(ev_update_finish);
1661 wallcycle_start(wcycle,ewcVSITECONSTR);
1664 shift_self(graph,state->box,state->x);
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);
1672 unshift_self(graph,state->box,state->x);
1674 wallcycle_stop(wcycle,ewcVSITECONSTR);
1677 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1678 if (ir->nstlist == -1 && bFirstIterate)
1680 gs.sig[eglsNABNSB] = nlh.nabnsb;
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,
1686 bFirstIterate ? &gs : NULL,
1687 (step_rel % gs.nstms == 0) &&
1688 (multisim_nsteps<0 || (step_rel<multisim_nsteps)),
1690 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
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)
1700 if (ir->nstlist == -1 && bFirstIterate)
1702 nlh.nabnsb = gs.set[eglsNABNSB];
1703 gs.set[eglsNABNSB] = 0;
1705 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1706 /* ############# END CALC EKIN AND PRESSURE ################# */
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. */
1714 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1715 trace(shake_vir),&tracevir))
1719 bFirstIterate = FALSE;
1722 /* only add constraint dvdl after constraints */
1723 enerd->term[F_DVDL_BONDED] += dvdl;
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);
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);
1732 /* ################# END UPDATE STEP 2 ################# */
1733 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1735 /* The coordinates (x) were unshifted in update */
1736 if (bFFscan && (shellfc==NULL || bConverged))
1738 if (print_forcefield(fplog,enerd->term,mdatoms->homenr,
1740 &(top_global->mols),mdatoms->massT,pres))
1744 fprintf(stderr,"\n");
1750 /* We will not sum ekinh_old,
1751 * so signal that we still have to do it.
1753 bSumEkinhOld = TRUE;
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
1764 /* Since this is called with the new coordinates state->x, I assume
1765 * we want the new box state->box too. / EL 20040121
1767 do_coupling(fplog,oenv,nfile,fnm,tcr,t,step,enerd->term,fr,
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);
1776 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1778 /* use the directly determined last velocity, not actually the averaged half steps */
1779 if (bTrotter && ir->eI==eiVV)
1781 enerd->term[F_EKIN] = last_ekin;
1783 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1787 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1791 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir,state,&MassQ);
1793 /* Check for excessively large energies */
1797 real etot_max = 1e200;
1799 real etot_max = 1e30;
1801 if (fabs(enerd->term[F_ETOT]) > etot_max)
1803 fprintf(stderr,"Energy too large (%g), giving up\n",
1804 enerd->term[F_ETOT]);
1807 /* ######### END PREPARING EDR OUTPUT ########### */
1809 /* Time for performance */
1810 if (((step % stepout) == 0) || bLastStep)
1812 runtime_upd_proc(runtime);
1818 gmx_bool do_dr,do_or;
1820 if (fplog && do_log && bDoExpanded)
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);
1826 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
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);
1838 upd_mdebin_step(mdebin);
1841 do_dr = do_per_step(step,ir->nstdisreout);
1842 do_or = do_per_step(step,ir->nstorireout);
1844 print_ebin(outf->fp_ene,do_ene,do_dr,do_or,do_log?fplog:NULL,
1846 eprNORMAL,bCompact,mdebin,fcd,groups,&(ir->opts));
1848 if (ir->ePull != epullNO)
1850 pull_print_output(ir->pull,step,t);
1853 if (do_per_step(step,ir->nstlog))
1855 if(fflush(fplog) != 0)
1857 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of disk space?");
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++)
1867 state->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
1870 /* Remaining runtime */
1871 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal() ))
1875 fprintf(stderr,"\n");
1877 print_time(stderr,runtime,step,ir,cr);
1880 /* Replica exchange */
1882 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1883 do_per_step(step,repl_ex_nst))
1885 bExchanged = replica_exchange(fplog,cr,repl_ex,
1889 if (bExchanged && DOMAINDECOMP(cr))
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,
1901 bStartingFromCpt = FALSE;
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.
1907 if ((state->flags & (1<<estPRES_PREV)) &&
1909 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1911 /* Store the pressure in t_state for pressure coupling
1912 * at the next MD step.
1914 copy_mat(pres,state->pres_prev);
1917 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1919 if ( (membed!=NULL) && (!bLastStep) )
1921 rescale_membed(step_rel,membed,state_global->x);
1928 /* read next frame from input trajectory */
1929 bNotLastFrame = read_next_frame(oenv,status,&rerun_fr);
1934 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
1938 if (!bRerunMD || !rerun_fr.bStep)
1940 /* increase the MD step number */
1945 cycles = wallcycle_stop(wcycle,ewcSTEP);
1946 if (DOMAINDECOMP(cr) && wcycle)
1948 dd_cycles_add(cr->dd,cycles,ddCyclStep);
1951 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1952 gs.set[eglsRESETCOUNTERS] != 0)
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;
1964 /* End of main MD loop */
1968 runtime_end(runtime);
1970 if (bRerunMD && MASTER(cr))
1975 if (!(cr->duty & DUTY_PME))
1977 /* Tell the PME only node to finish */
1983 if (ir->nstcalcenergy > 0 && !bRerunMD)
1985 print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
1986 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
1994 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
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);
2000 if (shellfc && fplog)
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",
2008 if (repl_ex_nst > 0 && MASTER(cr))
2010 print_replica_exchange_statistics(fplog,repl_ex);
2013 runtime->nsteps_done = step_rel;