2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
58 #include "md_support.h"
59 #include "md_logging.h"
74 #include "mpelogging.h"
76 #include "domdec_network.h"
82 #include "compute_io.h"
84 #include "checkpoint.h"
85 #include "mtop_util.h"
86 #include "sighandler.h"
89 #include "pme_loadbal.h"
92 #include "types/nlistheuristics.h"
93 #include "types/iteratedconstraints.h"
94 #include "nbnxn_cuda_data_mgmt.h"
104 #include "corewrap.h"
107 static void reset_all_counters(FILE *fplog,t_commrec *cr,
108 gmx_large_int_t step,
109 gmx_large_int_t *step_rel,t_inputrec *ir,
110 gmx_wallcycle_t wcycle,t_nrnb *nrnb,
111 gmx_runtime_t *runtime,
112 nbnxn_cuda_ptr_t cu_nbv)
114 char sbuf[STEPSTRSIZE];
116 /* Reset all the counters related to performance over the run */
117 md_print_warn(cr,fplog,"step %s: resetting all time and cycle counters\n",
118 gmx_step_str(step,sbuf));
122 nbnxn_cuda_reset_timings(cu_nbv);
125 wallcycle_stop(wcycle,ewcRUN);
126 wallcycle_reset_all(wcycle);
127 if (DOMAINDECOMP(cr))
129 reset_dd_statistics_counters(cr->dd);
132 ir->init_step += *step_rel;
133 ir->nsteps -= *step_rel;
135 wallcycle_start(wcycle,ewcRUN);
136 runtime_start(runtime);
137 print_date_and_time(fplog,cr->nodeid,"Restarted time",runtime);
140 double do_md(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
141 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
143 gmx_vsite_t *vsite,gmx_constr_t constr,
144 int stepout,t_inputrec *ir,
145 gmx_mtop_t *top_global,
147 t_state *state_global,
149 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
150 gmx_edsam_t ed,t_forcerec *fr,
151 int repl_ex_nst,int repl_ex_nex,int repl_ex_seed,gmx_membed_t membed,
152 real cpt_period,real max_hours,
153 const char *deviceOptions,
155 gmx_runtime_t *runtime)
158 gmx_large_int_t step,step_rel;
160 double t,t0,lam0[efptNR];
161 gmx_bool bGStatEveryStep,bGStat,bCalcVir,bCalcEner;
162 gmx_bool bNS,bNStList,bSimAnn,bStopCM,bRerunMD,bNotLastFrame=FALSE,
163 bFirstStep,bStateFromCP,bStateFromTPX,bInitStep,bLastStep,
164 bBornRadii,bStartingFromCpt;
165 gmx_bool bDoDHDL=FALSE,bDoFEP=FALSE,bDoExpanded=FALSE;
166 gmx_bool do_ene,do_log,do_verbose,bRerunWarnNoV=TRUE,
167 bForceUpdate=FALSE,bCPT;
169 gmx_bool bMasterState;
170 int force_flags,cglo_flags;
171 tensor force_vir,shake_vir,total_vir,tmp_vir,pres;
176 t_state *bufstate=NULL;
177 matrix *scale_tot,pcoupl_mu,M,ebox;
180 gmx_repl_ex_t repl_ex=NULL;
183 t_mdebin *mdebin=NULL;
184 df_history_t df_history;
189 gmx_enerdata_t *enerd;
191 gmx_global_stat_t gstat;
192 gmx_update_t upd=NULL;
195 gmx_rng_t mcrng=NULL;
197 gmx_groups_t *groups;
198 gmx_ekindata_t *ekind, *ekind_save;
199 gmx_shellfc_t shellfc;
200 int count,nconverged=0;
203 gmx_bool bIonize=FALSE;
204 gmx_bool bTCR=FALSE,bConverged=TRUE,bOK,bSumEkinhOld,bExchanged;
206 gmx_bool bResetCountersHalfMaxH=FALSE;
207 gmx_bool bVV,bIterativeCase,bFirstIterate,bTemp,bPres,bTrotter;
208 gmx_bool bUpdateDoLR;
211 atom_id *grpindex=NULL;
213 t_coupl_rec *tcr=NULL;
214 rvec *xcopy=NULL,*vcopy=NULL,*cbuf=NULL;
215 matrix boxcopy={{0}},lastbox;
217 real fom,oldfom,veta_save,pcurr,scalevir,tracevir;
224 real saved_conserved_quantity = 0;
229 char sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
230 int handled_stop_condition=gmx_stop_cond_none; /* compare to get_stop_condition*/
231 gmx_iterate_t iterate;
232 gmx_large_int_t multisim_nsteps=-1; /* number of steps to do before first multisim
233 simulation stops. If equal to zero, don't
234 communicate any more between multisims.*/
235 /* PME load balancing data for GPU kernels */
236 pme_load_balancing_t pme_loadbal=NULL;
238 gmx_bool bPMETuneTry=FALSE,bPMETuneRunning=FALSE;
241 /* Temporary addition for FAHCORE checkpointing */
245 /* Check for special mdrun options */
246 bRerunMD = (Flags & MD_RERUN);
247 bIonize = (Flags & MD_IONIZE);
248 bFFscan = (Flags & MD_FFSCAN);
249 bAppend = (Flags & MD_APPENDFILES);
250 if (Flags & MD_RESETCOUNTERSHALFWAY)
254 /* Signal to reset the counters half the simulation steps. */
255 wcycle_set_reset_counters(wcycle,ir->nsteps/2);
257 /* Signal to reset the counters halfway the simulation time. */
258 bResetCountersHalfMaxH = (max_hours > 0);
261 /* md-vv uses averaged full step velocities for T-control
262 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
263 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
265 if (bVV) /* to store the initial velocities while computing virial */
267 snew(cbuf,top_global->natoms);
269 /* all the iteratative cases - only if there are constraints */
270 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
271 gmx_iterate_init(&iterate,FALSE); /* The default value of iterate->bIterationActive is set to
272 false in this step. The correct value, true or false,
273 is set at each step, as it depends on the frequency of temperature
274 and pressure control.*/
275 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
279 /* Since we don't know if the frames read are related in any way,
280 * rebuild the neighborlist at every step.
283 ir->nstcalcenergy = 1;
287 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
289 nstglobalcomm = check_nstglobalcomm(fplog,cr,nstglobalcomm,ir);
290 bGStatEveryStep = (nstglobalcomm == 1);
292 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
295 "To reduce the energy communication with nstlist = -1\n"
296 "the neighbor list validity should not be checked at every step,\n"
297 "this means that exact integration is not guaranteed.\n"
298 "The neighbor list validity is checked after:\n"
299 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
300 "In most cases this will result in exact integration.\n"
301 "This reduces the energy communication by a factor of 2 to 3.\n"
302 "If you want less energy communication, set nstlist > 3.\n\n");
305 if (bRerunMD || bFFscan)
309 groups = &top_global->groups;
312 init_md(fplog,cr,ir,oenv,&t,&t0,state_global->lambda,
313 &(state_global->fep_state),lam0,
314 nrnb,top_global,&upd,
315 nfile,fnm,&outf,&mdebin,
316 force_vir,shake_vir,mu_tot,&bSimAnn,&vcm,state_global,Flags);
318 clear_mat(total_vir);
320 /* Energy terms and groups */
322 init_enerdata(top_global->groups.grps[egcENER].nr,ir->fepvals->n_lambda,
324 if (DOMAINDECOMP(cr))
330 snew(f,top_global->natoms);
333 /* lambda Monte carlo random number generator */
336 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
338 /* copy the state into df_history */
339 copy_df_history(&df_history,&state_global->dfhist);
341 /* Kinetic energy data */
343 init_ekindata(fplog,top_global,&(ir->opts),ekind);
344 /* needed for iteration of constraints */
346 init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
347 /* Copy the cos acceleration to the groups struct */
348 ekind->cosacc.cos_accel = ir->cos_accel;
350 gstat = global_stat_init(ir);
353 /* Check for polarizable models and flexible constraints */
354 shellfc = init_shell_flexcon(fplog,
355 top_global,n_flexible_constraints(constr),
356 (ir->bContinuation ||
357 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
358 NULL : state_global->x);
362 #ifdef GMX_THREAD_MPI
363 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
365 set_deform_reference_box(upd,
366 deform_init_init_step_tpx,
367 deform_init_box_tpx);
368 #ifdef GMX_THREAD_MPI
369 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
374 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
375 if ((io > 2000) && MASTER(cr))
377 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
381 if (DOMAINDECOMP(cr)) {
382 top = dd_init_local_top(top_global);
385 dd_init_local_state(cr->dd,state_global,state);
387 if (DDMASTER(cr->dd) && ir->nstfout) {
388 snew(f_global,state_global->natoms);
392 /* Initialize the particle decomposition and split the topology */
393 top = split_system(fplog,top_global,ir,cr);
395 pd_cg_range(cr,&fr->cg0,&fr->hcg);
396 pd_at_range(cr,&a0,&a1);
398 top = gmx_mtop_generate_local_top(top_global,ir);
401 a1 = top_global->natoms;
404 forcerec_set_excl_load(fr,top,cr);
406 state = partdec_init_local_state(cr,state_global);
409 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
412 set_vsite_top(vsite,top,mdatoms,cr);
415 if (ir->ePBC != epbcNONE && !fr->bMolPBC) {
416 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
420 make_local_shells(cr,mdatoms,shellfc);
423 init_bonded_thread_force_reduction(fr,&top->idef);
425 if (ir->pull && PAR(cr)) {
426 dd_make_local_pull_groups(NULL,ir->pull,mdatoms);
430 if (DOMAINDECOMP(cr))
432 /* Distribute the charge groups over the nodes from the master node */
433 dd_partition_system(fplog,ir->init_step,cr,TRUE,1,
434 state_global,top_global,ir,
435 state,&f,mdatoms,top,fr,
436 vsite,shellfc,constr,
441 update_mdatoms(mdatoms,state->lambda[efptMASS]);
443 if (opt2bSet("-cpi",nfile,fnm))
445 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr);
449 bStateFromCP = FALSE;
456 /* Update mdebin with energy history if appending to output files */
457 if ( Flags & MD_APPENDFILES )
459 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
463 /* We might have read an energy history from checkpoint,
464 * free the allocated memory and reset the counts.
466 done_energyhistory(&state_global->enerhist);
467 init_energyhistory(&state_global->enerhist);
470 /* Set the initial energy history in state by updating once */
471 update_energyhistory(&state_global->enerhist,mdebin);
474 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
476 /* Set the random state if we read a checkpoint file */
477 set_stochd_state(upd,state);
480 if (state->flags & (1<<estMC_RNG))
482 set_mc_state(mcrng,state);
485 /* Initialize constraints */
487 if (!DOMAINDECOMP(cr))
488 set_constraints(constr,top,ir,mdatoms,cr);
491 /* Check whether we have to GCT stuff */
492 bTCR = ftp2bSet(efGCT,nfile,fnm);
495 fprintf(stderr,"Will do General Coupling Theory!\n");
497 gnx = top_global->mols.nr;
499 for(i=0; (i<gnx); i++) {
506 /* We need to be sure replica exchange can only occur
507 * when the energies are current */
508 check_nst_param(fplog,cr,"nstcalcenergy",ir->nstcalcenergy,
509 "repl_ex_nst",&repl_ex_nst);
510 /* This check needs to happen before inter-simulation
511 * signals are initialized, too */
513 if (repl_ex_nst > 0 && MASTER(cr))
515 repl_ex = init_replica_exchange(fplog,cr->ms,state_global,ir,
516 repl_ex_nst,repl_ex_nex,repl_ex_seed);
519 /* PME tuning is only supported with GPUs or PME nodes and not with rerun */
520 if ((Flags & MD_TUNEPME) &&
521 EEL_PME(fr->eeltype) &&
522 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
525 pme_loadbal_init(&pme_loadbal,ir,state->box,fr->ic,fr->pmedata);
527 if (cr->duty & DUTY_PME)
529 /* Start tuning right away, as we can't measure the load */
530 bPMETuneRunning = TRUE;
534 /* Separate PME nodes, we can measure the PP/PME load balance */
539 if (!ir->bContinuation && !bRerunMD)
541 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
543 /* Set the velocities of frozen particles to zero */
544 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
548 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
558 /* Constrain the initial coordinates and velocities */
559 do_constrain_first(fplog,constr,ir,mdatoms,state,f,
560 graph,cr,nrnb,fr,top,shake_vir);
564 /* Construct the virtual sites for the initial configuration */
565 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
566 top->idef.iparams,top->idef.il,
567 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
573 /* set free energy calculation frequency as the minimum of nstdhdl, nstexpanded, and nstrepl_ex_nst*/
574 nstfep = ir->fepvals->nstdhdl;
575 if (ir->bExpanded && (nstfep > ir->expandedvals->nstexpanded))
577 nstfep = ir->expandedvals->nstexpanded;
579 if (repl_ex_nst > 0 && nstfep > repl_ex_nst)
581 nstfep = repl_ex_nst;
584 /* I'm assuming we need global communication the first time! MRS */
585 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
586 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM:0)
587 | (bVV ? CGLO_PRESSURE:0)
588 | (bVV ? CGLO_CONSTRAINT:0)
589 | (bRerunMD ? CGLO_RERUNMD:0)
590 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN:0));
592 bSumEkinhOld = FALSE;
593 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
594 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
595 constr,NULL,FALSE,state->box,
596 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,cglo_flags);
597 if (ir->eI == eiVVAK) {
598 /* a second call to get the half step temperature initialized as well */
599 /* we do the same call as above, but turn the pressure off -- internally to
600 compute_globals, this is recognized as a velocity verlet half-step
601 kinetic energy calculation. This minimized excess variables, but
602 perhaps loses some logic?*/
604 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
605 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
606 constr,NULL,FALSE,state->box,
607 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
608 cglo_flags &~ (CGLO_STOPCM | CGLO_PRESSURE));
611 /* Calculate the initial half step temperature, and save the ekinh_old */
612 if (!(Flags & MD_STARTFROMCPT))
614 for(i=0; (i<ir->opts.ngtc); i++)
616 copy_mat(ekind->tcstat[i].ekinh,ekind->tcstat[i].ekinh_old);
621 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
622 and there is no previous step */
625 /* if using an iterative algorithm, we need to create a working directory for the state. */
628 bufstate = init_bufstate(state);
632 snew(xcopy,state->natoms);
633 snew(vcopy,state->natoms);
634 copy_rvecn(state->x,xcopy,0,state->natoms);
635 copy_rvecn(state->v,vcopy,0,state->natoms);
636 copy_mat(state->box,boxcopy);
639 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
640 temperature control */
641 trotter_seq = init_npt_vars(ir,state,&MassQ,bTrotter);
645 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
648 "RMS relative constraint deviation after constraining: %.2e\n",
649 constr_rmsd(constr,FALSE));
651 if (EI_STATE_VELOCITY(ir->eI))
653 fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
657 fprintf(stderr,"starting md rerun '%s', reading coordinates from"
658 " input trajectory '%s'\n\n",
659 *(top_global->name),opt2fn("-rerun",nfile,fnm));
662 fprintf(stderr,"Calculated time to finish depends on nsteps from "
663 "run input file,\nwhich may not correspond to the time "
664 "needed to process input trajectory.\n\n");
670 fprintf(stderr,"starting mdrun '%s'\n",
671 *(top_global->name));
674 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
678 sprintf(tbuf,"%s","infinite");
680 if (ir->init_step > 0)
682 fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
683 gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
684 gmx_step_str(ir->init_step,sbuf2),
685 ir->init_step*ir->delta_t);
689 fprintf(stderr,"%s steps, %s ps.\n",
690 gmx_step_str(ir->nsteps,sbuf),tbuf);
696 /* Set and write start time */
697 runtime_start(runtime);
698 print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
699 wallcycle_start(wcycle,ewcRUN);
705 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
707 chkpt_ret=fcCheckPointParallel( cr->nodeid,
709 if ( chkpt_ret == 0 )
710 gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", 0 );
714 /***********************************************************
718 ************************************************************/
720 /* if rerunMD then read coordinates and velocities from input trajectory */
723 if (getenv("GMX_FORCE_UPDATE"))
731 bNotLastFrame = read_first_frame(oenv,&status,
732 opt2fn("-rerun",nfile,fnm),
733 &rerun_fr,TRX_NEED_X | TRX_READ_V);
734 if (rerun_fr.natoms != top_global->natoms)
737 "Number of atoms in trajectory (%d) does not match the "
738 "run input file (%d)\n",
739 rerun_fr.natoms,top_global->natoms);
741 if (ir->ePBC != epbcNONE)
745 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);
747 if (max_cutoff2(ir->ePBC,rerun_fr.box) < sqr(fr->rlistlong))
749 gmx_fatal(FARGS,"Rerun trajectory frame step %d time %f has too small box dimensions",rerun_fr.step,rerun_fr.time);
756 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
759 if (ir->ePBC != epbcNONE)
761 /* Set the shift vectors.
762 * Necessary here when have a static box different from the tpr box.
764 calc_shifts(rerun_fr.box,fr->shift_vec);
768 /* loop over MD steps or if rerunMD to end of input trajectory */
770 /* Skip the first Nose-Hoover integration when we get the state from tpx */
771 bStateFromTPX = !bStateFromCP;
772 bInitStep = bFirstStep && (bStateFromTPX || bVV);
773 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
775 bSumEkinhOld = FALSE;
778 init_global_signals(&gs,cr,ir,repl_ex_nst);
780 step = ir->init_step;
783 if (ir->nstlist == -1)
785 init_nlistheuristics(&nlh,bGStatEveryStep,step);
788 if (MULTISIM(cr) && (repl_ex_nst <=0 ))
790 /* check how many steps are left in other sims */
791 multisim_nsteps=get_multisim_nsteps(cr, ir->nsteps);
795 /* and stop now if we should */
796 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
797 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
798 while (!bLastStep || (bRerunMD && bNotLastFrame)) {
800 wallcycle_start(wcycle,ewcSTEP);
802 GMX_MPE_LOG(ev_timestep1);
805 if (rerun_fr.bStep) {
806 step = rerun_fr.step;
807 step_rel = step - ir->init_step;
809 if (rerun_fr.bTime) {
819 bLastStep = (step_rel == ir->nsteps);
820 t = t0 + step*ir->delta_t;
823 if (ir->efep != efepNO || ir->bSimTemp)
825 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
826 requiring different logic. */
828 set_current_lambdas(step,ir->fepvals,bRerunMD,&rerun_fr,state_global,state,lam0);
829 bDoDHDL = do_per_step(step,ir->fepvals->nstdhdl);
830 bDoFEP = (do_per_step(step,nstfep) && (ir->efep != efepNO));
831 bDoExpanded = (do_per_step(step,ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
836 update_annealing_target_temp(&(ir->opts),t);
841 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
843 for(i=0; i<state_global->natoms; i++)
845 copy_rvec(rerun_fr.x[i],state_global->x[i]);
849 for(i=0; i<state_global->natoms; i++)
851 copy_rvec(rerun_fr.v[i],state_global->v[i]);
856 for(i=0; i<state_global->natoms; i++)
858 clear_rvec(state_global->v[i]);
862 fprintf(stderr,"\nWARNING: Some frames do not contain velocities.\n"
863 " Ekin, temperature and pressure are incorrect,\n"
864 " the virial will be incorrect when constraints are present.\n"
866 bRerunWarnNoV = FALSE;
870 copy_mat(rerun_fr.box,state_global->box);
871 copy_mat(state_global->box,state->box);
873 if (vsite && (Flags & MD_RERUN_VSITE))
875 if (DOMAINDECOMP(cr))
877 gmx_fatal(FARGS,"Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
881 /* Following is necessary because the graph may get out of sync
882 * with the coordinates if we only have every N'th coordinate set
884 mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
885 shift_self(graph,state->box,state->x);
887 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
888 top->idef.iparams,top->idef.il,
889 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
892 unshift_self(graph,state->box,state->x);
897 /* Stop Center of Mass motion */
898 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step,ir->nstcomm));
900 /* Copy back starting coordinates in case we're doing a forcefield scan */
903 for(ii=0; (ii<state->natoms); ii++)
905 copy_rvec(xcopy[ii],state->x[ii]);
906 copy_rvec(vcopy[ii],state->v[ii]);
908 copy_mat(boxcopy,state->box);
913 /* for rerun MD always do Neighbour Searching */
914 bNS = (bFirstStep || ir->nstlist != 0);
919 /* Determine whether or not to do Neighbour Searching and LR */
920 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
922 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
923 (ir->nstlist == -1 && nlh.nabnsb > 0));
925 if (bNS && ir->nstlist == -1)
927 set_nlistheuristics(&nlh,bFirstStep || bExchanged || bDoFEP, step);
931 /* check whether we should stop because another simulation has
935 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
936 (multisim_nsteps != ir->nsteps) )
943 "Stopping simulation %d because another one has finished\n",
947 gs.sig[eglsCHKPT] = 1;
952 /* < 0 means stop at next step, > 0 means stop at next NS step */
953 if ( (gs.set[eglsSTOPCOND] < 0 ) ||
954 ( (gs.set[eglsSTOPCOND] > 0 ) && ( bNS || ir->nstlist==0)) )
959 /* Determine whether or not to update the Born radii if doing GB */
960 bBornRadii=bFirstStep;
961 if (ir->implicit_solvent && (step % ir->nstgbradii==0))
966 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
967 do_verbose = bVerbose &&
968 (step % stepout == 0 || bFirstStep || bLastStep);
970 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
978 bMasterState = FALSE;
979 /* Correct the new box if it is too skewed */
980 if (DYNAMIC_BOX(*ir))
982 if (correct_box(fplog,step,state->box,graph))
987 if (DOMAINDECOMP(cr) && bMasterState)
989 dd_collect_state(cr->dd,state,state_global);
993 if (DOMAINDECOMP(cr))
995 /* Repartition the domain decomposition */
996 wallcycle_start(wcycle,ewcDOMDEC);
997 dd_partition_system(fplog,step,cr,
998 bMasterState,nstglobalcomm,
999 state_global,top_global,ir,
1000 state,&f,mdatoms,top,fr,
1001 vsite,shellfc,constr,
1003 do_verbose && !bPMETuneRunning);
1004 wallcycle_stop(wcycle,ewcDOMDEC);
1005 /* If using an iterative integrator, reallocate space to match the decomposition */
1009 if (MASTER(cr) && do_log && !bFFscan)
1011 print_ebin_header(fplog,step,t,state->lambda[efptFEP]); /* can we improve the information printed here? */
1014 if (ir->efep != efepNO)
1016 update_mdatoms(mdatoms,state->lambda[efptMASS]);
1019 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1022 /* We need the kinetic energy at minus the half step for determining
1023 * the full step kinetic energy and possibly for T-coupling.*/
1024 /* This may not be quite working correctly yet . . . . */
1025 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1026 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1027 constr,NULL,FALSE,state->box,
1028 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1029 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1031 clear_mat(force_vir);
1033 /* Ionize the atoms if necessary */
1036 ionize(fplog,oenv,mdatoms,top_global,t,ir,state->x,state->v,
1037 mdatoms->start,mdatoms->start+mdatoms->homenr,state->box,cr);
1040 /* Update force field in ffscan program */
1043 if (update_forcefield(fplog,
1045 mdatoms->nr,state->x,state->box))
1053 GMX_MPE_LOG(ev_timestep2);
1055 /* We write a checkpoint at this MD step when:
1056 * either at an NS step when we signalled through gs,
1057 * or at the last step (but not when we do not want confout),
1058 * but never at the first step or with rerun.
1060 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1061 (bLastStep && (Flags & MD_CONFOUT))) &&
1062 step > ir->init_step && !bRerunMD);
1065 gs.set[eglsCHKPT] = 0;
1068 /* Determine the energy and pressure:
1069 * at nstcalcenergy steps and at energy output steps (set below).
1071 if (EI_VV(ir->eI) && (!bInitStep))
1073 /* for vv, the first half of the integration actually corresponds
1074 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1075 but the virial needs to be calculated on both the current step and the 'next' step. Future
1076 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1078 bCalcEner = do_per_step(step-1,ir->nstcalcenergy);
1079 bCalcVir = bCalcEner ||
1080 (ir->epc != epcNO && (do_per_step(step,ir->nstpcouple) || do_per_step(step-1,ir->nstpcouple)));
1084 bCalcEner = do_per_step(step,ir->nstcalcenergy);
1085 bCalcVir = bCalcEner ||
1086 (ir->epc != epcNO && do_per_step(step,ir->nstpcouple));
1089 /* Do we need global communication ? */
1090 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1091 do_per_step(step,nstglobalcomm) ||
1092 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1094 do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
1096 if (do_ene || do_log)
1103 /* these CGLO_ options remain the same throughout the iteration */
1104 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1105 (bGStat ? CGLO_GSTAT : 0)
1108 force_flags = (GMX_FORCE_STATECHANGED |
1109 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1110 GMX_FORCE_ALLFORCES |
1112 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1113 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1114 (bDoFEP ? GMX_FORCE_DHDL : 0)
1119 if(do_per_step(step,ir->nstcalclr))
1121 force_flags |= GMX_FORCE_DO_LR;
1127 /* Now is the time to relax the shells */
1128 count=relax_shell_flexcon(fplog,cr,bVerbose,bFFscan ? step+1 : step,
1130 bStopCM,top,top_global,
1132 state,f,force_vir,mdatoms,
1133 nrnb,wcycle,graph,groups,
1134 shellfc,fr,bBornRadii,t,mu_tot,
1135 state->natoms,&bConverged,vsite,
1146 /* The coordinates (x) are shifted (to get whole molecules)
1148 * This is parallellized as well, and does communication too.
1149 * Check comments in sim_util.c
1151 do_force(fplog,cr,ir,step,nrnb,wcycle,top,top_global,groups,
1152 state->box,state->x,&state->hist,
1153 f,force_vir,mdatoms,enerd,fcd,
1154 state->lambda,graph,
1155 fr,vsite,mu_tot,t,outf->fp_field,ed,bBornRadii,
1156 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1159 GMX_BARRIER(cr->mpi_comm_mygroup);
1163 mu_aver = calc_mu_aver(cr,state->x,mdatoms->chargeA,
1164 mu_tot,&top_global->mols,mdatoms,gnx,grpindex);
1167 if (bTCR && bFirstStep)
1169 tcr=init_coupling(fplog,nfile,fnm,cr,fr,mdatoms,&(top->idef));
1170 fprintf(fplog,"Done init_coupling\n");
1174 if (bVV && !bStartingFromCpt && !bRerunMD)
1175 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1177 if (ir->eI==eiVV && bInitStep)
1179 /* if using velocity verlet with full time step Ekin,
1180 * take the first half step only to compute the
1181 * virial for the first step. From there,
1182 * revert back to the initial coordinates
1183 * so that the input is actually the initial step.
1185 copy_rvecn(state->v,cbuf,0,state->natoms); /* should make this better for parallelizing? */
1187 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1188 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ1);
1191 /* If we are using twin-range interactions where the long-range component
1192 * is only evaluated every nstcalclr>1 steps, we should do a special update
1193 * step to combine the long-range forces on these steps.
1194 * For nstcalclr=1 this is not done, since the forces would have been added
1195 * directly to the short-range forces already.
1197 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1199 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,
1200 f,bUpdateDoLR,fr->f_twin,fcd,
1201 ekind,M,wcycle,upd,bInitStep,etrtVELOCITY1,
1202 cr,nrnb,constr,&top->idef);
1204 if (bIterativeCase && do_per_step(step-1,ir->nstpcouple) && !bInitStep)
1206 gmx_iterate_init(&iterate,TRUE);
1208 /* for iterations, we save these vectors, as we will be self-consistently iterating
1211 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1213 /* save the state */
1214 if (iterate.bIterationActive) {
1215 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1218 bFirstIterate = TRUE;
1219 while (bFirstIterate || iterate.bIterationActive)
1221 if (iterate.bIterationActive)
1223 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1224 if (bFirstIterate && bTrotter)
1226 /* The first time through, we need a decent first estimate
1227 of veta(t+dt) to compute the constraints. Do
1228 this by computing the box volume part of the
1229 trotter integration at this time. Nothing else
1230 should be changed by this routine here. If
1231 !(first time), we start with the previous value
1234 veta_save = state->veta;
1235 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ0);
1236 vetanew = state->veta;
1237 state->veta = veta_save;
1242 if ( !bRerunMD || rerun_fr.bV || bForceUpdate) { /* Why is rerun_fr.bV here? Unclear. */
1245 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1246 state,fr->bMolPBC,graph,f,
1247 &top->idef,shake_vir,NULL,
1248 cr,nrnb,wcycle,upd,constr,
1249 bInitStep,TRUE,bCalcVir,vetanew);
1251 if (!bOK && !bFFscan)
1253 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1258 { /* Need to unshift here if a do_force has been
1259 called in the previous step */
1260 unshift_self(graph,state->box,state->x);
1263 /* if VV, compute the pressure and constraints */
1264 /* For VV2, we strictly only need this if using pressure
1265 * control, but we really would like to have accurate pressures
1267 * Think about ways around this in the future?
1268 * For now, keep this choice in comments.
1270 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1271 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1273 bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK));
1274 if (bCalcEner && ir->eI==eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1276 bSumEkinhOld = TRUE;
1278 /* for vv, the first half of the integration actually corresponds to the previous step.
1279 So we need information from the last step in the first half of the integration */
1280 if (bGStat || do_per_step(step-1,nstglobalcomm)) {
1281 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1282 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1283 constr,NULL,FALSE,state->box,
1284 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1287 | (bTemp ? CGLO_TEMPERATURE:0)
1288 | (bPres ? CGLO_PRESSURE : 0)
1289 | (bPres ? CGLO_CONSTRAINT : 0)
1290 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1291 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1294 /* explanation of above:
1295 a) We compute Ekin at the full time step
1296 if 1) we are using the AveVel Ekin, and it's not the
1297 initial step, or 2) if we are using AveEkin, but need the full
1298 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1299 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1300 EkinAveVel because it's needed for the pressure */
1302 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1307 m_add(force_vir,shake_vir,total_vir); /* we need the un-dispersion corrected total vir here */
1308 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ2);
1315 /* We need the kinetic energy at minus the half step for determining
1316 * the full step kinetic energy and possibly for T-coupling.*/
1317 /* This may not be quite working correctly yet . . . . */
1318 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1319 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1320 constr,NULL,FALSE,state->box,
1321 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1322 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1327 if (iterate.bIterationActive &&
1328 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1329 state->veta,&vetanew))
1333 bFirstIterate = FALSE;
1336 if (bTrotter && !bInitStep) {
1337 enerd->term[F_DVDL_BONDED] += dvdl; /* only add after iterations */
1338 copy_mat(shake_vir,state->svir_prev);
1339 copy_mat(force_vir,state->fvir_prev);
1340 if (IR_NVT_TROTTER(ir) && ir->eI==eiVV) {
1341 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1342 enerd->term[F_TEMP] = sum_ekin(&(ir->opts),ekind,NULL,(ir->eI==eiVV),FALSE,FALSE);
1343 enerd->term[F_EKIN] = trace(ekind->ekin);
1346 /* if it's the initial step, we performed this first step just to get the constraint virial */
1347 if (bInitStep && ir->eI==eiVV) {
1348 copy_rvecn(cbuf,state->v,0,state->natoms);
1351 if (fr->bSepDVDL && fplog && do_log)
1353 fprintf(fplog,sepdvdlformat,"Constraint",0.0,dvdl);
1355 enerd->term[F_DVDL_BONDED] += dvdl;
1357 GMX_MPE_LOG(ev_timestep1);
1360 /* MRS -- now done iterating -- compute the conserved quantity */
1362 saved_conserved_quantity = compute_conserved_from_auxiliary(ir,state,&MassQ);
1365 last_ekin = enerd->term[F_EKIN];
1367 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1369 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1371 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1374 sum_dhdl(enerd,state->lambda,ir->fepvals);
1378 /* ######## END FIRST UPDATE STEP ############## */
1379 /* ######## If doing VV, we now have v(dt) ###### */
1381 /* perform extended ensemble sampling in lambda - we don't
1382 actually move to the new state before outputting
1383 statistics, but if performing simulated tempering, we
1384 do update the velocities and the tau_t. */
1386 lamnew = ExpandedEnsembleDynamics(fplog,ir,enerd,state,&MassQ,&df_history,step,mcrng,state->v,mdatoms);
1388 /* ################## START TRAJECTORY OUTPUT ################# */
1390 /* Now we have the energies and forces corresponding to the
1391 * coordinates at time t. We must output all of this before
1393 * for RerunMD t is read from input trajectory
1395 GMX_MPE_LOG(ev_output_start);
1398 if (do_per_step(step,ir->nstxout)) { mdof_flags |= MDOF_X; }
1399 if (do_per_step(step,ir->nstvout)) { mdof_flags |= MDOF_V; }
1400 if (do_per_step(step,ir->nstfout)) { mdof_flags |= MDOF_F; }
1401 if (do_per_step(step,ir->nstxtcout)) { mdof_flags |= MDOF_XTC; }
1402 if (bCPT) { mdof_flags |= MDOF_CPT; };
1404 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1407 /* Enforce writing positions and velocities at end of run */
1408 mdof_flags |= (MDOF_X | MDOF_V);
1413 fcReportProgress( ir->nsteps, step );
1415 /* sync bCPT and fc record-keeping */
1416 if (bCPT && MASTER(cr))
1417 fcRequestCheckPoint();
1420 if (mdof_flags != 0)
1422 wallcycle_start(wcycle,ewcTRAJ);
1425 if (state->flags & (1<<estLD_RNG))
1427 get_stochd_state(upd,state);
1429 if (state->flags & (1<<estMC_RNG))
1431 get_mc_state(mcrng,state);
1437 state_global->ekinstate.bUpToDate = FALSE;
1441 update_ekinstate(&state_global->ekinstate,ekind);
1442 state_global->ekinstate.bUpToDate = TRUE;
1444 update_energyhistory(&state_global->enerhist,mdebin);
1445 if (ir->efep!=efepNO || ir->bSimTemp)
1447 state_global->fep_state = state->fep_state; /* MRS: seems kludgy. The code should be
1448 structured so this isn't necessary.
1449 Note this reassignment is only necessary
1450 for single threads.*/
1451 copy_df_history(&state_global->dfhist,&df_history);
1455 write_traj(fplog,cr,outf,mdof_flags,top_global,
1456 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
1463 if (bLastStep && step_rel == ir->nsteps &&
1464 (Flags & MD_CONFOUT) && MASTER(cr) &&
1465 !bRerunMD && !bFFscan)
1467 /* x and v have been collected in write_traj,
1468 * because a checkpoint file will always be written
1471 fprintf(stderr,"\nWriting final coordinates.\n");
1474 /* Make molecules whole only for confout writing */
1475 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
1477 write_sto_conf_mtop(ftp2fn(efSTO,nfile,fnm),
1478 *top_global->name,top_global,
1479 state_global->x,state_global->v,
1480 ir->ePBC,state->box);
1483 wallcycle_stop(wcycle,ewcTRAJ);
1485 GMX_MPE_LOG(ev_output_finish);
1487 /* kludge -- virial is lost with restart for NPT control. Must restart */
1488 if (bStartingFromCpt && bVV)
1490 copy_mat(state->svir_prev,shake_vir);
1491 copy_mat(state->fvir_prev,force_vir);
1493 /* ################## END TRAJECTORY OUTPUT ################ */
1495 /* Determine the wallclock run time up till now */
1496 run_time = gmx_gettime() - (double)runtime->real;
1498 /* Check whether everything is still allright */
1499 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1500 #ifdef GMX_THREAD_MPI
1505 /* this is just make gs.sig compatible with the hack
1506 of sending signals around by MPI_Reduce with together with
1508 if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
1509 gs.sig[eglsSTOPCOND]=1;
1510 if ( gmx_get_stop_condition() == gmx_stop_cond_next )
1511 gs.sig[eglsSTOPCOND]=-1;
1512 /* < 0 means stop at next step, > 0 means stop at next NS step */
1516 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1517 gmx_get_signal_name(),
1518 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1522 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1523 gmx_get_signal_name(),
1524 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1526 handled_stop_condition=(int)gmx_get_stop_condition();
1528 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1529 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1530 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1532 /* Signal to terminate the run */
1533 gs.sig[eglsSTOPCOND] = 1;
1536 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1538 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1541 if (bResetCountersHalfMaxH && MASTER(cr) &&
1542 run_time > max_hours*60.0*60.0*0.495)
1544 gs.sig[eglsRESETCOUNTERS] = 1;
1547 if (ir->nstlist == -1 && !bRerunMD)
1549 /* When bGStatEveryStep=FALSE, global_stat is only called
1550 * when we check the atom displacements, not at NS steps.
1551 * This means that also the bonded interaction count check is not
1552 * performed immediately after NS. Therefore a few MD steps could
1553 * be performed with missing interactions.
1554 * But wrong energies are never written to file,
1555 * since energies are only written after global_stat
1558 if (step >= nlh.step_nscheck)
1560 nlh.nabnsb = natoms_beyond_ns_buffer(ir,fr,&top->cgs,
1561 nlh.scale_tot,state->x);
1565 /* This is not necessarily true,
1566 * but step_nscheck is determined quite conservatively.
1572 /* In parallel we only have to check for checkpointing in steps
1573 * where we do global communication,
1574 * otherwise the other nodes don't know.
1576 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1579 run_time >= nchkpt*cpt_period*60.0)) &&
1580 gs.set[eglsCHKPT] == 0)
1582 gs.sig[eglsCHKPT] = 1;
1585 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1590 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1592 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1594 gmx_bool bIfRandomize;
1595 bIfRandomize = update_randomize_velocities(ir,step,mdatoms,state,upd,&top->idef,constr);
1596 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1597 if (constr && bIfRandomize)
1599 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1600 state,fr->bMolPBC,graph,f,
1601 &top->idef,tmp_vir,NULL,
1602 cr,nrnb,wcycle,upd,constr,
1603 bInitStep,TRUE,bCalcVir,vetanew);
1608 if (bIterativeCase && do_per_step(step,ir->nstpcouple))
1610 gmx_iterate_init(&iterate,TRUE);
1611 /* for iterations, we save these vectors, as we will be redoing the calculations */
1612 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1615 bFirstIterate = TRUE;
1616 while (bFirstIterate || iterate.bIterationActive)
1618 /* We now restore these vectors to redo the calculation with improved extended variables */
1619 if (iterate.bIterationActive)
1621 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1624 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1625 so scroll down for that logic */
1627 /* ######### START SECOND UPDATE STEP ################# */
1628 GMX_MPE_LOG(ev_update_start);
1629 /* Box is changed in update() when we do pressure coupling,
1630 * but we should still use the old box for energy corrections and when
1631 * writing it to the energy file, so it matches the trajectory files for
1632 * the same timestep above. Make a copy in a separate array.
1634 copy_mat(state->box,lastbox);
1637 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1639 wallcycle_start(wcycle,ewcUPDATE);
1641 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1644 if (iterate.bIterationActive)
1652 /* we use a new value of scalevir to converge the iterations faster */
1653 scalevir = tracevir/trace(shake_vir);
1655 msmul(shake_vir,scalevir,shake_vir);
1656 m_add(force_vir,shake_vir,total_vir);
1657 clear_mat(shake_vir);
1659 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ3);
1660 /* We can only do Berendsen coupling after we have summed
1661 * the kinetic energy or virial. Since the happens
1662 * in global_state after update, we should only do it at
1663 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1668 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1669 update_pcouple(fplog,step,ir,state,pcoupl_mu,M,wcycle,
1675 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1677 /* velocity half-step update */
1678 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1679 bUpdateDoLR,fr->f_twin,fcd,
1680 ekind,M,wcycle,upd,FALSE,etrtVELOCITY2,
1681 cr,nrnb,constr,&top->idef);
1684 /* Above, initialize just copies ekinh into ekin,
1685 * it doesn't copy position (for VV),
1686 * and entire integrator for MD.
1691 copy_rvecn(state->x,cbuf,0,state->natoms);
1693 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1695 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1696 bUpdateDoLR,fr->f_twin,fcd,
1697 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1698 wallcycle_stop(wcycle,ewcUPDATE);
1700 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,state,
1701 fr->bMolPBC,graph,f,
1702 &top->idef,shake_vir,force_vir,
1703 cr,nrnb,wcycle,upd,constr,
1704 bInitStep,FALSE,bCalcVir,state->veta);
1708 /* erase F_EKIN and F_TEMP here? */
1709 /* just compute the kinetic energy at the half step to perform a trotter step */
1710 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1711 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1712 constr,NULL,FALSE,lastbox,
1713 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1714 cglo_flags | CGLO_TEMPERATURE
1716 wallcycle_start(wcycle,ewcUPDATE);
1717 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ4);
1718 /* now we know the scaling, we can compute the positions again again */
1719 copy_rvecn(cbuf,state->x,0,state->natoms);
1721 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1723 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1724 bUpdateDoLR,fr->f_twin,fcd,
1725 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1726 wallcycle_stop(wcycle,ewcUPDATE);
1728 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1729 /* are the small terms in the shake_vir here due
1730 * to numerical errors, or are they important
1731 * physically? I'm thinking they are just errors, but not completely sure.
1732 * For now, will call without actually constraining, constr=NULL*/
1733 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1734 state,fr->bMolPBC,graph,f,
1735 &top->idef,tmp_vir,force_vir,
1736 cr,nrnb,wcycle,upd,NULL,
1737 bInitStep,FALSE,bCalcVir,
1740 if (!bOK && !bFFscan)
1742 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1745 if (fr->bSepDVDL && fplog && do_log)
1747 fprintf(fplog,sepdvdlformat,"Constraint dV/dl",0.0,dvdl);
1749 enerd->term[F_DVDL_BONDED] += dvdl;
1753 /* Need to unshift here */
1754 unshift_self(graph,state->box,state->x);
1757 GMX_BARRIER(cr->mpi_comm_mygroup);
1758 GMX_MPE_LOG(ev_update_finish);
1762 wallcycle_start(wcycle,ewcVSITECONSTR);
1765 shift_self(graph,state->box,state->x);
1767 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
1768 top->idef.iparams,top->idef.il,
1769 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
1773 unshift_self(graph,state->box,state->x);
1775 wallcycle_stop(wcycle,ewcVSITECONSTR);
1778 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1779 /* With Leap-Frog we can skip compute_globals at
1780 * non-communication steps, but we need to calculate
1781 * the kinetic energy one step before communication.
1783 if (bGStat || (!EI_VV(ir->eI)&&do_per_step(step+1,nstglobalcomm)))
1785 if (ir->nstlist == -1 && bFirstIterate)
1787 gs.sig[eglsNABNSB] = nlh.nabnsb;
1789 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1790 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1792 bFirstIterate ? &gs : NULL,
1793 (step_rel % gs.nstms == 0) &&
1794 (multisim_nsteps<0 || (step_rel<multisim_nsteps)),
1796 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1798 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1799 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1800 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1801 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1802 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1803 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1806 if (ir->nstlist == -1 && bFirstIterate)
1808 nlh.nabnsb = gs.set[eglsNABNSB];
1809 gs.set[eglsNABNSB] = 0;
1812 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1813 /* ############# END CALC EKIN AND PRESSURE ################# */
1815 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1816 the virial that should probably be addressed eventually. state->veta has better properies,
1817 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1818 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1820 if (iterate.bIterationActive &&
1821 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1822 trace(shake_vir),&tracevir))
1826 bFirstIterate = FALSE;
1829 /* only add constraint dvdl after constraints */
1830 enerd->term[F_DVDL_BONDED] += dvdl;
1831 if (!bVV || bRerunMD)
1833 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1834 sum_dhdl(enerd,state->lambda,ir->fepvals);
1836 update_box(fplog,step,ir,mdatoms,state,graph,f,
1837 ir->nstlist==-1 ? &nlh.scale_tot : NULL,pcoupl_mu,nrnb,wcycle,upd,bInitStep,FALSE);
1839 /* ################# END UPDATE STEP 2 ################# */
1840 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1842 /* The coordinates (x) were unshifted in update */
1843 if (bFFscan && (shellfc==NULL || bConverged))
1845 if (print_forcefield(fplog,enerd->term,mdatoms->homenr,
1847 &(top_global->mols),mdatoms->massT,pres))
1851 fprintf(stderr,"\n");
1857 /* We will not sum ekinh_old,
1858 * so signal that we still have to do it.
1860 bSumEkinhOld = TRUE;
1865 /* Only do GCT when the relaxation of shells (minimization) has converged,
1866 * otherwise we might be coupling to bogus energies.
1867 * In parallel we must always do this, because the other sims might
1871 /* Since this is called with the new coordinates state->x, I assume
1872 * we want the new box state->box too. / EL 20040121
1874 do_coupling(fplog,oenv,nfile,fnm,tcr,t,step,enerd->term,fr,
1876 mdatoms,&(top->idef),mu_aver,
1877 top_global->mols.nr,cr,
1878 state->box,total_vir,pres,
1879 mu_tot,state->x,f,bConverged);
1883 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1885 /* use the directly determined last velocity, not actually the averaged half steps */
1886 if (bTrotter && ir->eI==eiVV)
1888 enerd->term[F_EKIN] = last_ekin;
1890 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1894 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1898 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir,state,&MassQ);
1900 /* Check for excessively large energies */
1904 real etot_max = 1e200;
1906 real etot_max = 1e30;
1908 if (fabs(enerd->term[F_ETOT]) > etot_max)
1910 fprintf(stderr,"Energy too large (%g), giving up\n",
1911 enerd->term[F_ETOT]);
1914 /* ######### END PREPARING EDR OUTPUT ########### */
1916 /* Time for performance */
1917 if (((step % stepout) == 0) || bLastStep)
1919 runtime_upd_proc(runtime);
1925 gmx_bool do_dr,do_or;
1927 if (fplog && do_log && bDoExpanded)
1929 /* only needed if doing expanded ensemble */
1930 PrintFreeEnergyInfoToFile(fplog,ir->fepvals,ir->expandedvals,ir->bSimTemp?ir->simtempvals:NULL,
1931 &df_history,state->fep_state,ir->nstlog,step);
1933 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1937 upd_mdebin(mdebin,bDoDHDL, TRUE,
1938 t,mdatoms->tmass,enerd,state,
1939 ir->fepvals,ir->expandedvals,lastbox,
1940 shake_vir,force_vir,total_vir,pres,
1941 ekind,mu_tot,constr);
1945 upd_mdebin_step(mdebin);
1948 do_dr = do_per_step(step,ir->nstdisreout);
1949 do_or = do_per_step(step,ir->nstorireout);
1951 print_ebin(outf->fp_ene,do_ene,do_dr,do_or,do_log?fplog:NULL,
1953 eprNORMAL,bCompact,mdebin,fcd,groups,&(ir->opts));
1955 if (ir->ePull != epullNO)
1957 pull_print_output(ir->pull,step,t);
1960 if (do_per_step(step,ir->nstlog))
1962 if(fflush(fplog) != 0)
1964 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of disk space?");
1970 /* Have to do this part after outputting the logfile and the edr file */
1971 state->fep_state = lamnew;
1972 for (i=0;i<efptNR;i++)
1974 state_global->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
1977 /* Remaining runtime */
1978 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1982 fprintf(stderr,"\n");
1984 print_time(stderr,runtime,step,ir,cr);
1987 /* Replica exchange */
1989 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1990 do_per_step(step,repl_ex_nst))
1992 bExchanged = replica_exchange(fplog,cr,repl_ex,
1996 if (bExchanged && DOMAINDECOMP(cr))
1998 dd_partition_system(fplog,step,cr,TRUE,1,
1999 state_global,top_global,ir,
2000 state,&f,mdatoms,top,fr,
2001 vsite,shellfc,constr,
2008 bStartingFromCpt = FALSE;
2010 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
2011 /* With all integrators, except VV, we need to retain the pressure
2012 * at the current step for coupling at the next step.
2014 if ((state->flags & (1<<estPRES_PREV)) &&
2016 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
2018 /* Store the pressure in t_state for pressure coupling
2019 * at the next MD step.
2021 copy_mat(pres,state->pres_prev);
2024 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2026 if ( (membed!=NULL) && (!bLastStep) )
2028 rescale_membed(step_rel,membed,state_global->x);
2035 /* read next frame from input trajectory */
2036 bNotLastFrame = read_next_frame(oenv,status,&rerun_fr);
2041 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
2045 if (!bRerunMD || !rerun_fr.bStep)
2047 /* increase the MD step number */
2052 cycles = wallcycle_stop(wcycle,ewcSTEP);
2053 if (DOMAINDECOMP(cr) && wcycle)
2055 dd_cycles_add(cr->dd,cycles,ddCyclStep);
2058 if (bPMETuneRunning || bPMETuneTry)
2060 /* PME grid + cut-off optimization with GPUs or PME nodes */
2062 /* Count the total cycles over the last steps */
2063 cycles_pmes += cycles;
2065 /* We can only switch cut-off at NS steps */
2066 if (step % ir->nstlist == 0)
2068 /* PME grid + cut-off optimization with GPUs or PME nodes */
2071 if (DDMASTER(cr->dd))
2073 /* PME node load is too high, start tuning */
2074 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2076 dd_bcast(cr->dd,sizeof(gmx_bool),&bPMETuneRunning);
2078 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2080 bPMETuneTry = FALSE;
2083 if (bPMETuneRunning)
2085 /* init_step might not be a multiple of nstlist,
2086 * but the first cycle is always skipped anyhow.
2089 pme_load_balance(pme_loadbal,cr,
2090 (bVerbose && MASTER(cr)) ? stderr : NULL,
2092 ir,state,cycles_pmes,
2093 fr->ic,fr->nbv,&fr->pmedata,
2096 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2097 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2098 fr->rlist = fr->ic->rlist;
2099 fr->rlistlong = fr->ic->rlistlong;
2100 fr->rcoulomb = fr->ic->rcoulomb;
2101 fr->rvdw = fr->ic->rvdw;
2107 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2108 gs.set[eglsRESETCOUNTERS] != 0)
2110 /* Reset all the counters related to performance over the run */
2111 reset_all_counters(fplog,cr,step,&step_rel,ir,wcycle,nrnb,runtime,
2112 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2113 wcycle_set_reset_counters(wcycle,-1);
2114 /* Correct max_hours for the elapsed time */
2115 max_hours -= run_time/(60.0*60.0);
2116 bResetCountersHalfMaxH = FALSE;
2117 gs.set[eglsRESETCOUNTERS] = 0;
2121 /* End of main MD loop */
2125 runtime_end(runtime);
2127 if (bRerunMD && MASTER(cr))
2132 if (!(cr->duty & DUTY_PME))
2134 /* Tell the PME only node to finish */
2135 gmx_pme_send_finish(cr);
2140 if (ir->nstcalcenergy > 0 && !bRerunMD)
2142 print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
2143 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
2151 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2153 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)));
2154 fprintf(fplog,"Average number of atoms that crossed the half buffer length: %.1f\n\n",nlh.ab/nlh.nns);
2157 if (pme_loadbal != NULL)
2159 pme_loadbal_done(pme_loadbal,fplog);
2162 if (shellfc && fplog)
2164 fprintf(fplog,"Fraction of iterations that converged: %.2f %%\n",
2165 (nconverged*100.0)/step_rel);
2166 fprintf(fplog,"Average number of force evaluations per MD step: %.2f\n\n",
2170 if (repl_ex_nst > 0 && MASTER(cr))
2172 print_replica_exchange_statistics(fplog,repl_ex);
2175 runtime->nsteps_done = step_rel;