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, 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,bIterations,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;
225 real saved_conserved_quantity = 0;
230 char sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
231 int handled_stop_condition=gmx_stop_cond_none; /* compare to get_stop_condition*/
232 gmx_iterate_t iterate;
233 gmx_large_int_t multisim_nsteps=-1; /* number of steps to do before first multisim
234 simulation stops. If equal to zero, don't
235 communicate any more between multisims.*/
236 /* PME load balancing data for GPU kernels */
237 pme_load_balancing_t pme_loadbal=NULL;
239 gmx_bool bPMETuneTry=FALSE,bPMETuneRunning=FALSE;
242 /* Temporary addition for FAHCORE checkpointing */
246 /* Check for special mdrun options */
247 bRerunMD = (Flags & MD_RERUN);
248 bIonize = (Flags & MD_IONIZE);
249 bFFscan = (Flags & MD_FFSCAN);
250 bAppend = (Flags & MD_APPENDFILES);
251 if (Flags & MD_RESETCOUNTERSHALFWAY)
255 /* Signal to reset the counters half the simulation steps. */
256 wcycle_set_reset_counters(wcycle,ir->nsteps/2);
258 /* Signal to reset the counters halfway the simulation time. */
259 bResetCountersHalfMaxH = (max_hours > 0);
262 /* md-vv uses averaged full step velocities for T-control
263 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
264 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
266 if (bVV) /* to store the initial velocities while computing virial */
268 snew(cbuf,top_global->natoms);
270 /* all the iteratative cases - only if there are constraints */
271 bIterations = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
272 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
276 /* Since we don't know if the frames read are related in any way,
277 * rebuild the neighborlist at every step.
280 ir->nstcalcenergy = 1;
284 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
286 nstglobalcomm = check_nstglobalcomm(fplog,cr,nstglobalcomm,ir);
287 bGStatEveryStep = (nstglobalcomm == 1);
289 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
292 "To reduce the energy communication with nstlist = -1\n"
293 "the neighbor list validity should not be checked at every step,\n"
294 "this means that exact integration is not guaranteed.\n"
295 "The neighbor list validity is checked after:\n"
296 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
297 "In most cases this will result in exact integration.\n"
298 "This reduces the energy communication by a factor of 2 to 3.\n"
299 "If you want less energy communication, set nstlist > 3.\n\n");
302 if (bRerunMD || bFFscan)
306 groups = &top_global->groups;
309 init_md(fplog,cr,ir,oenv,&t,&t0,state_global->lambda,
310 &(state_global->fep_state),lam0,
311 nrnb,top_global,&upd,
312 nfile,fnm,&outf,&mdebin,
313 force_vir,shake_vir,mu_tot,&bSimAnn,&vcm,state_global,Flags);
315 clear_mat(total_vir);
317 /* Energy terms and groups */
319 init_enerdata(top_global->groups.grps[egcENER].nr,ir->fepvals->n_lambda,
321 if (DOMAINDECOMP(cr))
327 snew(f,top_global->natoms);
330 /* lambda Monte carlo random number generator */
333 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
335 /* copy the state into df_history */
336 copy_df_history(&df_history,&state_global->dfhist);
338 /* Kinetic energy data */
340 init_ekindata(fplog,top_global,&(ir->opts),ekind);
341 /* needed for iteration of constraints */
343 init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
344 /* Copy the cos acceleration to the groups struct */
345 ekind->cosacc.cos_accel = ir->cos_accel;
347 gstat = global_stat_init(ir);
350 /* Check for polarizable models and flexible constraints */
351 shellfc = init_shell_flexcon(fplog,
352 top_global,n_flexible_constraints(constr),
353 (ir->bContinuation ||
354 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
355 NULL : state_global->x);
359 #ifdef GMX_THREAD_MPI
360 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
362 set_deform_reference_box(upd,
363 deform_init_init_step_tpx,
364 deform_init_box_tpx);
365 #ifdef GMX_THREAD_MPI
366 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
371 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
372 if ((io > 2000) && MASTER(cr))
374 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
378 if (DOMAINDECOMP(cr)) {
379 top = dd_init_local_top(top_global);
382 dd_init_local_state(cr->dd,state_global,state);
384 if (DDMASTER(cr->dd) && ir->nstfout) {
385 snew(f_global,state_global->natoms);
389 /* Initialize the particle decomposition and split the topology */
390 top = split_system(fplog,top_global,ir,cr);
392 pd_cg_range(cr,&fr->cg0,&fr->hcg);
393 pd_at_range(cr,&a0,&a1);
395 top = gmx_mtop_generate_local_top(top_global,ir);
398 a1 = top_global->natoms;
401 forcerec_set_excl_load(fr,top,cr);
403 state = partdec_init_local_state(cr,state_global);
406 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
409 set_vsite_top(vsite,top,mdatoms,cr);
412 if (ir->ePBC != epbcNONE && !fr->bMolPBC) {
413 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
417 make_local_shells(cr,mdatoms,shellfc);
420 init_bonded_thread_force_reduction(fr,&top->idef);
422 if (ir->pull && PAR(cr)) {
423 dd_make_local_pull_groups(NULL,ir->pull,mdatoms);
427 if (DOMAINDECOMP(cr))
429 /* Distribute the charge groups over the nodes from the master node */
430 dd_partition_system(fplog,ir->init_step,cr,TRUE,1,
431 state_global,top_global,ir,
432 state,&f,mdatoms,top,fr,
433 vsite,shellfc,constr,
438 update_mdatoms(mdatoms,state->lambda[efptMASS]);
440 if (opt2bSet("-cpi",nfile,fnm))
442 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr);
446 bStateFromCP = FALSE;
453 /* Update mdebin with energy history if appending to output files */
454 if ( Flags & MD_APPENDFILES )
456 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
460 /* We might have read an energy history from checkpoint,
461 * free the allocated memory and reset the counts.
463 done_energyhistory(&state_global->enerhist);
464 init_energyhistory(&state_global->enerhist);
467 /* Set the initial energy history in state by updating once */
468 update_energyhistory(&state_global->enerhist,mdebin);
471 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
473 /* Set the random state if we read a checkpoint file */
474 set_stochd_state(upd,state);
477 if (state->flags & (1<<estMC_RNG))
479 set_mc_state(mcrng,state);
482 /* Initialize constraints */
484 if (!DOMAINDECOMP(cr))
485 set_constraints(constr,top,ir,mdatoms,cr);
488 /* Check whether we have to GCT stuff */
489 bTCR = ftp2bSet(efGCT,nfile,fnm);
492 fprintf(stderr,"Will do General Coupling Theory!\n");
494 gnx = top_global->mols.nr;
496 for(i=0; (i<gnx); i++) {
503 /* We need to be sure replica exchange can only occur
504 * when the energies are current */
505 check_nst_param(fplog,cr,"nstcalcenergy",ir->nstcalcenergy,
506 "repl_ex_nst",&repl_ex_nst);
507 /* This check needs to happen before inter-simulation
508 * signals are initialized, too */
510 if (repl_ex_nst > 0 && MASTER(cr))
512 repl_ex = init_replica_exchange(fplog,cr->ms,state_global,ir,
513 repl_ex_nst,repl_ex_nex,repl_ex_seed);
516 /* PME tuning is only supported with GPUs or PME nodes and not with rerun */
517 if ((Flags & MD_TUNEPME) &&
518 EEL_PME(fr->eeltype) &&
519 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
522 pme_loadbal_init(&pme_loadbal,ir,state->box,fr->ic,fr->pmedata);
524 if (cr->duty & DUTY_PME)
526 /* Start tuning right away, as we can't measure the load */
527 bPMETuneRunning = TRUE;
531 /* Separate PME nodes, we can measure the PP/PME load balance */
536 if (!ir->bContinuation && !bRerunMD)
538 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
540 /* Set the velocities of frozen particles to zero */
541 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
545 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
555 /* Constrain the initial coordinates and velocities */
556 do_constrain_first(fplog,constr,ir,mdatoms,state,f,
557 graph,cr,nrnb,fr,top,shake_vir);
561 /* Construct the virtual sites for the initial configuration */
562 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
563 top->idef.iparams,top->idef.il,
564 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
570 /* set free energy calculation frequency as the minimum of nstdhdl, nstexpanded, and nstrepl_ex_nst*/
571 nstfep = ir->fepvals->nstdhdl;
572 if (ir->bExpanded && (nstfep > ir->expandedvals->nstexpanded))
574 nstfep = ir->expandedvals->nstexpanded;
576 if (repl_ex_nst > 0 && repl_ex_nst > nstfep)
578 nstfep = repl_ex_nst;
581 /* I'm assuming we need global communication the first time! MRS */
582 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
583 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM:0)
584 | (bVV ? CGLO_PRESSURE:0)
585 | (bVV ? CGLO_CONSTRAINT:0)
586 | (bRerunMD ? CGLO_RERUNMD:0)
587 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN:0));
589 bSumEkinhOld = FALSE;
590 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
591 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
592 constr,NULL,FALSE,state->box,
593 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,cglo_flags);
594 if (ir->eI == eiVVAK) {
595 /* a second call to get the half step temperature initialized as well */
596 /* we do the same call as above, but turn the pressure off -- internally to
597 compute_globals, this is recognized as a velocity verlet half-step
598 kinetic energy calculation. This minimized excess variables, but
599 perhaps loses some logic?*/
601 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
602 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
603 constr,NULL,FALSE,state->box,
604 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
605 cglo_flags &~ (CGLO_STOPCM | CGLO_PRESSURE));
608 /* Calculate the initial half step temperature, and save the ekinh_old */
609 if (!(Flags & MD_STARTFROMCPT))
611 for(i=0; (i<ir->opts.ngtc); i++)
613 copy_mat(ekind->tcstat[i].ekinh,ekind->tcstat[i].ekinh_old);
618 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
619 and there is no previous step */
622 /* if using an iterative algorithm, we need to create a working directory for the state. */
625 bufstate = init_bufstate(state);
629 snew(xcopy,state->natoms);
630 snew(vcopy,state->natoms);
631 copy_rvecn(state->x,xcopy,0,state->natoms);
632 copy_rvecn(state->v,vcopy,0,state->natoms);
633 copy_mat(state->box,boxcopy);
636 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
637 temperature control */
638 trotter_seq = init_npt_vars(ir,state,&MassQ,bTrotter);
642 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
645 "RMS relative constraint deviation after constraining: %.2e\n",
646 constr_rmsd(constr,FALSE));
648 if (EI_STATE_VELOCITY(ir->eI))
650 fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
654 fprintf(stderr,"starting md rerun '%s', reading coordinates from"
655 " input trajectory '%s'\n\n",
656 *(top_global->name),opt2fn("-rerun",nfile,fnm));
659 fprintf(stderr,"Calculated time to finish depends on nsteps from "
660 "run input file,\nwhich may not correspond to the time "
661 "needed to process input trajectory.\n\n");
667 fprintf(stderr,"starting mdrun '%s'\n",
668 *(top_global->name));
671 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
675 sprintf(tbuf,"%s","infinite");
677 if (ir->init_step > 0)
679 fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
680 gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
681 gmx_step_str(ir->init_step,sbuf2),
682 ir->init_step*ir->delta_t);
686 fprintf(stderr,"%s steps, %s ps.\n",
687 gmx_step_str(ir->nsteps,sbuf),tbuf);
693 /* Set and write start time */
694 runtime_start(runtime);
695 print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
696 wallcycle_start(wcycle,ewcRUN);
702 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
704 chkpt_ret=fcCheckPointParallel( cr->nodeid,
706 if ( chkpt_ret == 0 )
707 gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", 0 );
711 /***********************************************************
715 ************************************************************/
717 /* if rerunMD then read coordinates and velocities from input trajectory */
720 if (getenv("GMX_FORCE_UPDATE"))
728 bNotLastFrame = read_first_frame(oenv,&status,
729 opt2fn("-rerun",nfile,fnm),
730 &rerun_fr,TRX_NEED_X | TRX_READ_V);
731 if (rerun_fr.natoms != top_global->natoms)
734 "Number of atoms in trajectory (%d) does not match the "
735 "run input file (%d)\n",
736 rerun_fr.natoms,top_global->natoms);
738 if (ir->ePBC != epbcNONE)
742 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);
744 if (max_cutoff2(ir->ePBC,rerun_fr.box) < sqr(fr->rlistlong))
746 gmx_fatal(FARGS,"Rerun trajectory frame step %d time %f has too small box dimensions",rerun_fr.step,rerun_fr.time);
753 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
756 if (ir->ePBC != epbcNONE)
758 /* Set the shift vectors.
759 * Necessary here when have a static box different from the tpr box.
761 calc_shifts(rerun_fr.box,fr->shift_vec);
765 /* loop over MD steps or if rerunMD to end of input trajectory */
767 /* Skip the first Nose-Hoover integration when we get the state from tpx */
768 bStateFromTPX = !bStateFromCP;
769 bInitStep = bFirstStep && (bStateFromTPX || bVV);
770 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
772 bSumEkinhOld = FALSE;
775 init_global_signals(&gs,cr,ir,repl_ex_nst);
777 step = ir->init_step;
780 if (ir->nstlist == -1)
782 init_nlistheuristics(&nlh,bGStatEveryStep,step);
785 if (MULTISIM(cr) && (repl_ex_nst <=0 ))
787 /* check how many steps are left in other sims */
788 multisim_nsteps=get_multisim_nsteps(cr, ir->nsteps);
792 /* and stop now if we should */
793 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
794 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
795 while (!bLastStep || (bRerunMD && bNotLastFrame)) {
797 wallcycle_start(wcycle,ewcSTEP);
799 GMX_MPE_LOG(ev_timestep1);
802 if (rerun_fr.bStep) {
803 step = rerun_fr.step;
804 step_rel = step - ir->init_step;
806 if (rerun_fr.bTime) {
816 bLastStep = (step_rel == ir->nsteps);
817 t = t0 + step*ir->delta_t;
820 if (ir->efep != efepNO || ir->bSimTemp)
822 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
823 requiring different logic. */
825 set_current_lambdas(step,ir->fepvals,bRerunMD,&rerun_fr,state_global,state,lam0);
826 bDoDHDL = do_per_step(step,ir->fepvals->nstdhdl);
827 bDoFEP = (do_per_step(step,nstfep) && (ir->efep != efepNO));
828 bDoExpanded = (do_per_step(step,ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
833 update_annealing_target_temp(&(ir->opts),t);
838 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
840 for(i=0; i<state_global->natoms; i++)
842 copy_rvec(rerun_fr.x[i],state_global->x[i]);
846 for(i=0; i<state_global->natoms; i++)
848 copy_rvec(rerun_fr.v[i],state_global->v[i]);
853 for(i=0; i<state_global->natoms; i++)
855 clear_rvec(state_global->v[i]);
859 fprintf(stderr,"\nWARNING: Some frames do not contain velocities.\n"
860 " Ekin, temperature and pressure are incorrect,\n"
861 " the virial will be incorrect when constraints are present.\n"
863 bRerunWarnNoV = FALSE;
867 copy_mat(rerun_fr.box,state_global->box);
868 copy_mat(state_global->box,state->box);
870 if (vsite && (Flags & MD_RERUN_VSITE))
872 if (DOMAINDECOMP(cr))
874 gmx_fatal(FARGS,"Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
878 /* Following is necessary because the graph may get out of sync
879 * with the coordinates if we only have every N'th coordinate set
881 mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
882 shift_self(graph,state->box,state->x);
884 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
885 top->idef.iparams,top->idef.il,
886 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
889 unshift_self(graph,state->box,state->x);
894 /* Stop Center of Mass motion */
895 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step,ir->nstcomm));
897 /* Copy back starting coordinates in case we're doing a forcefield scan */
900 for(ii=0; (ii<state->natoms); ii++)
902 copy_rvec(xcopy[ii],state->x[ii]);
903 copy_rvec(vcopy[ii],state->v[ii]);
905 copy_mat(boxcopy,state->box);
910 /* for rerun MD always do Neighbour Searching */
911 bNS = (bFirstStep || ir->nstlist != 0);
916 /* Determine whether or not to do Neighbour Searching and LR */
917 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
919 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
920 (ir->nstlist == -1 && nlh.nabnsb > 0));
922 if (bNS && ir->nstlist == -1)
924 set_nlistheuristics(&nlh,bFirstStep || bExchanged || bDoFEP, step);
928 /* check whether we should stop because another simulation has
932 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
933 (multisim_nsteps != ir->nsteps) )
940 "Stopping simulation %d because another one has finished\n",
944 gs.sig[eglsCHKPT] = 1;
949 /* < 0 means stop at next step, > 0 means stop at next NS step */
950 if ( (gs.set[eglsSTOPCOND] < 0 ) ||
951 ( (gs.set[eglsSTOPCOND] > 0 ) && ( bNS || ir->nstlist==0)) )
956 /* Determine whether or not to update the Born radii if doing GB */
957 bBornRadii=bFirstStep;
958 if (ir->implicit_solvent && (step % ir->nstgbradii==0))
963 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
964 do_verbose = bVerbose &&
965 (step % stepout == 0 || bFirstStep || bLastStep);
967 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
975 bMasterState = FALSE;
976 /* Correct the new box if it is too skewed */
977 if (DYNAMIC_BOX(*ir))
979 if (correct_box(fplog,step,state->box,graph))
984 if (DOMAINDECOMP(cr) && bMasterState)
986 dd_collect_state(cr->dd,state,state_global);
990 if (DOMAINDECOMP(cr))
992 /* Repartition the domain decomposition */
993 wallcycle_start(wcycle,ewcDOMDEC);
994 dd_partition_system(fplog,step,cr,
995 bMasterState,nstglobalcomm,
996 state_global,top_global,ir,
997 state,&f,mdatoms,top,fr,
998 vsite,shellfc,constr,
1000 do_verbose && !bPMETuneRunning);
1001 wallcycle_stop(wcycle,ewcDOMDEC);
1002 /* If using an iterative integrator, reallocate space to match the decomposition */
1006 if (MASTER(cr) && do_log && !bFFscan)
1008 print_ebin_header(fplog,step,t,state->lambda[efptFEP]); /* can we improve the information printed here? */
1011 if (ir->efep != efepNO)
1013 update_mdatoms(mdatoms,state->lambda[efptMASS]);
1016 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1019 /* We need the kinetic energy at minus the half step for determining
1020 * the full step kinetic energy and possibly for T-coupling.*/
1021 /* This may not be quite working correctly yet . . . . */
1022 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1023 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1024 constr,NULL,FALSE,state->box,
1025 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1026 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1028 clear_mat(force_vir);
1030 /* Ionize the atoms if necessary */
1033 ionize(fplog,oenv,mdatoms,top_global,t,ir,state->x,state->v,
1034 mdatoms->start,mdatoms->start+mdatoms->homenr,state->box,cr);
1037 /* Update force field in ffscan program */
1040 if (update_forcefield(fplog,
1042 mdatoms->nr,state->x,state->box))
1050 GMX_MPE_LOG(ev_timestep2);
1052 /* We write a checkpoint at this MD step when:
1053 * either at an NS step when we signalled through gs,
1054 * or at the last step (but not when we do not want confout),
1055 * but never at the first step or with rerun.
1057 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1058 (bLastStep && (Flags & MD_CONFOUT))) &&
1059 step > ir->init_step && !bRerunMD);
1062 gs.set[eglsCHKPT] = 0;
1065 /* Determine the energy and pressure:
1066 * at nstcalcenergy steps and at energy output steps (set below).
1068 if (EI_VV(ir->eI) && (!bInitStep))
1070 /* for vv, the first half actually corresponds to the last step */
1071 bCalcEner = do_per_step(step-1,ir->nstcalcenergy);
1075 bCalcEner = do_per_step(step,ir->nstcalcenergy);
1077 bCalcVir = bCalcEner ||
1078 (ir->epc != epcNO && do_per_step(step,ir->nstpcouple));
1080 /* Do we need global communication ? */
1081 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1082 do_per_step(step,nstglobalcomm) ||
1083 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1085 do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
1087 if (do_ene || do_log)
1094 /* these CGLO_ options remain the same throughout the iteration */
1095 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1096 (bGStat ? CGLO_GSTAT : 0)
1099 force_flags = (GMX_FORCE_STATECHANGED |
1100 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1101 GMX_FORCE_ALLFORCES |
1103 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1104 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1105 (bDoFEP ? GMX_FORCE_DHDL : 0)
1110 if(do_per_step(step,ir->nstcalclr))
1112 force_flags |= GMX_FORCE_DO_LR;
1118 /* Now is the time to relax the shells */
1119 count=relax_shell_flexcon(fplog,cr,bVerbose,bFFscan ? step+1 : step,
1121 bStopCM,top,top_global,
1123 state,f,force_vir,mdatoms,
1124 nrnb,wcycle,graph,groups,
1125 shellfc,fr,bBornRadii,t,mu_tot,
1126 state->natoms,&bConverged,vsite,
1137 /* The coordinates (x) are shifted (to get whole molecules)
1139 * This is parallellized as well, and does communication too.
1140 * Check comments in sim_util.c
1142 do_force(fplog,cr,ir,step,nrnb,wcycle,top,top_global,groups,
1143 state->box,state->x,&state->hist,
1144 f,force_vir,mdatoms,enerd,fcd,
1145 state->lambda,graph,
1146 fr,vsite,mu_tot,t,outf->fp_field,ed,bBornRadii,
1147 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1150 GMX_BARRIER(cr->mpi_comm_mygroup);
1154 mu_aver = calc_mu_aver(cr,state->x,mdatoms->chargeA,
1155 mu_tot,&top_global->mols,mdatoms,gnx,grpindex);
1158 if (bTCR && bFirstStep)
1160 tcr=init_coupling(fplog,nfile,fnm,cr,fr,mdatoms,&(top->idef));
1161 fprintf(fplog,"Done init_coupling\n");
1165 if (bVV && !bStartingFromCpt && !bRerunMD)
1166 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1168 if (ir->eI==eiVV && bInitStep)
1170 /* if using velocity verlet with full time step Ekin,
1171 * take the first half step only to compute the
1172 * virial for the first step. From there,
1173 * revert back to the initial coordinates
1174 * so that the input is actually the initial step.
1176 copy_rvecn(state->v,cbuf,0,state->natoms); /* should make this better for parallelizing? */
1178 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1179 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ1);
1182 /* If we are using twin-range interactions where the long-range component
1183 * is only evaluated every nstcalclr>1 steps, we should do a special update
1184 * step to combine the long-range forces on these steps.
1185 * For nstcalclr=1 this is not done, since the forces would have been added
1186 * directly to the short-range forces already.
1188 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1190 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,
1191 f,bUpdateDoLR,fr->f_twin,fcd,
1192 ekind,M,wcycle,upd,bInitStep,etrtVELOCITY1,
1193 cr,nrnb,constr,&top->idef);
1197 gmx_iterate_init(&iterate,bIterations && !bInitStep);
1199 /* for iterations, we save these vectors, as we will be self-consistently iterating
1202 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1204 /* save the state */
1205 if (bIterations && iterate.bIterate) {
1206 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1209 bFirstIterate = TRUE;
1210 while (bFirstIterate || (bIterations && iterate.bIterate))
1212 if (bIterations && iterate.bIterate)
1214 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1215 if (bFirstIterate && bTrotter)
1217 /* The first time through, we need a decent first estimate
1218 of veta(t+dt) to compute the constraints. Do
1219 this by computing the box volume part of the
1220 trotter integration at this time. Nothing else
1221 should be changed by this routine here. If
1222 !(first time), we start with the previous value
1225 veta_save = state->veta;
1226 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ0);
1227 vetanew = state->veta;
1228 state->veta = veta_save;
1233 if ( !bRerunMD || rerun_fr.bV || bForceUpdate) { /* Why is rerun_fr.bV here? Unclear. */
1236 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1237 state,fr->bMolPBC,graph,f,
1238 &top->idef,shake_vir,NULL,
1239 cr,nrnb,wcycle,upd,constr,
1240 bInitStep,TRUE,bCalcVir,vetanew);
1242 if (!bOK && !bFFscan)
1244 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1249 { /* Need to unshift here if a do_force has been
1250 called in the previous step */
1251 unshift_self(graph,state->box,state->x);
1255 /* if VV, compute the pressure and constraints */
1256 /* For VV2, we strictly only need this if using pressure
1257 * control, but we really would like to have accurate pressures
1259 * Think about ways around this in the future?
1260 * For now, keep this choice in comments.
1262 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1263 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1265 bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK));
1266 if (bCalcEner && ir->eI==eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1268 bSumEkinhOld = TRUE;
1270 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1271 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1272 constr,NULL,FALSE,state->box,
1273 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1276 | (bStopCM ? CGLO_STOPCM : 0)
1277 | (bTemp ? CGLO_TEMPERATURE:0)
1278 | (bPres ? CGLO_PRESSURE : 0)
1279 | (bPres ? CGLO_CONSTRAINT : 0)
1280 | ((bIterations && iterate.bIterate) ? CGLO_ITERATE : 0)
1281 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1284 /* explanation of above:
1285 a) We compute Ekin at the full time step
1286 if 1) we are using the AveVel Ekin, and it's not the
1287 initial step, or 2) if we are using AveEkin, but need the full
1288 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1289 b) If we are using EkinAveEkin for the kinetic energy for the temperture control, we still feed in
1290 EkinAveVel because it's needed for the pressure */
1292 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1297 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ2);
1304 /* We need the kinetic energy at minus the half step for determining
1305 * the full step kinetic energy and possibly for T-coupling.*/
1306 /* This may not be quite working correctly yet . . . . */
1307 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1308 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1309 constr,NULL,FALSE,state->box,
1310 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1311 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1315 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1320 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1321 state->veta,&vetanew))
1325 bFirstIterate = FALSE;
1328 if (bTrotter && !bInitStep) {
1329 enerd->term[F_DVDL_BONDED] += dvdl; /* only add after iterations */
1330 copy_mat(shake_vir,state->svir_prev);
1331 copy_mat(force_vir,state->fvir_prev);
1332 if (IR_NVT_TROTTER(ir) && ir->eI==eiVV) {
1333 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1334 enerd->term[F_TEMP] = sum_ekin(&(ir->opts),ekind,NULL,(ir->eI==eiVV),FALSE,FALSE);
1335 enerd->term[F_EKIN] = trace(ekind->ekin);
1338 /* if it's the initial step, we performed this first step just to get the constraint virial */
1339 if (bInitStep && ir->eI==eiVV) {
1340 copy_rvecn(cbuf,state->v,0,state->natoms);
1343 if (fr->bSepDVDL && fplog && do_log)
1345 fprintf(fplog,sepdvdlformat,"Constraint",0.0,dvdl);
1347 enerd->term[F_DVDL_BONDED] += dvdl;
1349 GMX_MPE_LOG(ev_timestep1);
1352 /* MRS -- now done iterating -- compute the conserved quantity */
1354 saved_conserved_quantity = compute_conserved_from_auxiliary(ir,state,&MassQ);
1357 last_ekin = enerd->term[F_EKIN];
1359 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1361 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1363 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1364 sum_dhdl(enerd,state->lambda,ir->fepvals);
1367 /* ######## END FIRST UPDATE STEP ############## */
1368 /* ######## If doing VV, we now have v(dt) ###### */
1370 /* perform extended ensemble sampling in lambda - we don't
1371 actually move to the new state before outputting
1372 statistics, but if performing simulated tempering, we
1373 do update the velocities and the tau_t. */
1375 lamnew = ExpandedEnsembleDynamics(fplog,ir,enerd,state,&MassQ,&df_history,step,mcrng,state->v,mdatoms);
1377 /* ################## START TRAJECTORY OUTPUT ################# */
1379 /* Now we have the energies and forces corresponding to the
1380 * coordinates at time t. We must output all of this before
1382 * for RerunMD t is read from input trajectory
1384 GMX_MPE_LOG(ev_output_start);
1387 if (do_per_step(step,ir->nstxout)) { mdof_flags |= MDOF_X; }
1388 if (do_per_step(step,ir->nstvout)) { mdof_flags |= MDOF_V; }
1389 if (do_per_step(step,ir->nstfout)) { mdof_flags |= MDOF_F; }
1390 if (do_per_step(step,ir->nstxtcout)) { mdof_flags |= MDOF_XTC; }
1391 if (bCPT) { mdof_flags |= MDOF_CPT; };
1393 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1396 /* Enforce writing positions and velocities at end of run */
1397 mdof_flags |= (MDOF_X | MDOF_V);
1402 fcReportProgress( ir->nsteps, step );
1404 /* sync bCPT and fc record-keeping */
1405 if (bCPT && MASTER(cr))
1406 fcRequestCheckPoint();
1409 if (mdof_flags != 0)
1411 wallcycle_start(wcycle,ewcTRAJ);
1414 if (state->flags & (1<<estLD_RNG))
1416 get_stochd_state(upd,state);
1418 if (state->flags & (1<<estMC_RNG))
1420 get_mc_state(mcrng,state);
1426 state_global->ekinstate.bUpToDate = FALSE;
1430 update_ekinstate(&state_global->ekinstate,ekind);
1431 state_global->ekinstate.bUpToDate = TRUE;
1433 update_energyhistory(&state_global->enerhist,mdebin);
1434 if (ir->efep!=efepNO || ir->bSimTemp)
1436 state_global->fep_state = state->fep_state; /* MRS: seems kludgy. The code should be
1437 structured so this isn't necessary.
1438 Note this reassignment is only necessary
1439 for single threads.*/
1440 copy_df_history(&state_global->dfhist,&df_history);
1444 write_traj(fplog,cr,outf,mdof_flags,top_global,
1445 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
1452 if (bLastStep && step_rel == ir->nsteps &&
1453 (Flags & MD_CONFOUT) && MASTER(cr) &&
1454 !bRerunMD && !bFFscan)
1456 /* x and v have been collected in write_traj,
1457 * because a checkpoint file will always be written
1460 fprintf(stderr,"\nWriting final coordinates.\n");
1463 /* Make molecules whole only for confout writing */
1464 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
1466 write_sto_conf_mtop(ftp2fn(efSTO,nfile,fnm),
1467 *top_global->name,top_global,
1468 state_global->x,state_global->v,
1469 ir->ePBC,state->box);
1472 wallcycle_stop(wcycle,ewcTRAJ);
1474 GMX_MPE_LOG(ev_output_finish);
1476 /* kludge -- virial is lost with restart for NPT control. Must restart */
1477 if (bStartingFromCpt && bVV)
1479 copy_mat(state->svir_prev,shake_vir);
1480 copy_mat(state->fvir_prev,force_vir);
1482 /* ################## END TRAJECTORY OUTPUT ################ */
1484 /* Determine the wallclock run time up till now */
1485 run_time = gmx_gettime() - (double)runtime->real;
1487 /* Check whether everything is still allright */
1488 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1489 #ifdef GMX_THREAD_MPI
1494 /* this is just make gs.sig compatible with the hack
1495 of sending signals around by MPI_Reduce with together with
1497 if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
1498 gs.sig[eglsSTOPCOND]=1;
1499 if ( gmx_get_stop_condition() == gmx_stop_cond_next )
1500 gs.sig[eglsSTOPCOND]=-1;
1501 /* < 0 means stop at next step, > 0 means stop at next NS step */
1505 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1506 gmx_get_signal_name(),
1507 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1511 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1512 gmx_get_signal_name(),
1513 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1515 handled_stop_condition=(int)gmx_get_stop_condition();
1517 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1518 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1519 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1521 /* Signal to terminate the run */
1522 gs.sig[eglsSTOPCOND] = 1;
1525 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1527 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1530 if (bResetCountersHalfMaxH && MASTER(cr) &&
1531 run_time > max_hours*60.0*60.0*0.495)
1533 gs.sig[eglsRESETCOUNTERS] = 1;
1536 if (ir->nstlist == -1 && !bRerunMD)
1538 /* When bGStatEveryStep=FALSE, global_stat is only called
1539 * when we check the atom displacements, not at NS steps.
1540 * This means that also the bonded interaction count check is not
1541 * performed immediately after NS. Therefore a few MD steps could
1542 * be performed with missing interactions.
1543 * But wrong energies are never written to file,
1544 * since energies are only written after global_stat
1547 if (step >= nlh.step_nscheck)
1549 nlh.nabnsb = natoms_beyond_ns_buffer(ir,fr,&top->cgs,
1550 nlh.scale_tot,state->x);
1554 /* This is not necessarily true,
1555 * but step_nscheck is determined quite conservatively.
1561 /* In parallel we only have to check for checkpointing in steps
1562 * where we do global communication,
1563 * otherwise the other nodes don't know.
1565 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1568 run_time >= nchkpt*cpt_period*60.0)) &&
1569 gs.set[eglsCHKPT] == 0)
1571 gs.sig[eglsCHKPT] = 1;
1575 /* at the start of step, randomize the velocities */
1576 if (ETC_ANDERSEN(ir->etc) && EI_VV(ir->eI))
1578 gmx_bool bDoAndersenConstr;
1579 bDoAndersenConstr = (constr && update_randomize_velocities(ir,step,mdatoms,state,upd,&top->idef,constr));
1580 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1581 if (bDoAndersenConstr)
1583 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1584 state,fr->bMolPBC,graph,f,
1585 &top->idef,tmp_vir,NULL,
1586 cr,nrnb,wcycle,upd,constr,
1587 bInitStep,TRUE,bCalcVir,vetanew);
1593 gmx_iterate_init(&iterate,bIterations);
1596 /* for iterations, we save these vectors, as we will be redoing the calculations */
1597 if (bIterations && iterate.bIterate)
1599 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1601 bFirstIterate = TRUE;
1602 while (bFirstIterate || (bIterations && iterate.bIterate))
1604 /* We now restore these vectors to redo the calculation with improved extended variables */
1607 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1610 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1611 so scroll down for that logic */
1613 /* ######### START SECOND UPDATE STEP ################# */
1614 GMX_MPE_LOG(ev_update_start);
1615 /* Box is changed in update() when we do pressure coupling,
1616 * but we should still use the old box for energy corrections and when
1617 * writing it to the energy file, so it matches the trajectory files for
1618 * the same timestep above. Make a copy in a separate array.
1620 copy_mat(state->box,lastbox);
1623 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1625 wallcycle_start(wcycle,ewcUPDATE);
1627 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1630 if (bIterations && iterate.bIterate)
1638 /* we use a new value of scalevir to converge the iterations faster */
1639 scalevir = tracevir/trace(shake_vir);
1641 msmul(shake_vir,scalevir,shake_vir);
1642 m_add(force_vir,shake_vir,total_vir);
1643 clear_mat(shake_vir);
1645 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ3);
1646 /* We can only do Berendsen coupling after we have summed
1647 * the kinetic energy or virial. Since the happens
1648 * in global_state after update, we should only do it at
1649 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1654 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1655 update_pcouple(fplog,step,ir,state,pcoupl_mu,M,wcycle,
1661 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1663 /* velocity half-step update */
1664 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1665 bUpdateDoLR,fr->f_twin,fcd,
1666 ekind,M,wcycle,upd,FALSE,etrtVELOCITY2,
1667 cr,nrnb,constr,&top->idef);
1670 /* Above, initialize just copies ekinh into ekin,
1671 * it doesn't copy position (for VV),
1672 * and entire integrator for MD.
1677 copy_rvecn(state->x,cbuf,0,state->natoms);
1679 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1681 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1682 bUpdateDoLR,fr->f_twin,fcd,
1683 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1684 wallcycle_stop(wcycle,ewcUPDATE);
1686 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,state,
1687 fr->bMolPBC,graph,f,
1688 &top->idef,shake_vir,force_vir,
1689 cr,nrnb,wcycle,upd,constr,
1690 bInitStep,FALSE,bCalcVir,state->veta);
1694 /* erase F_EKIN and F_TEMP here? */
1695 /* just compute the kinetic energy at the half step to perform a trotter step */
1696 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1697 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1698 constr,NULL,FALSE,lastbox,
1699 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1700 cglo_flags | CGLO_TEMPERATURE
1702 wallcycle_start(wcycle,ewcUPDATE);
1703 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ4);
1704 /* now we know the scaling, we can compute the positions again again */
1705 copy_rvecn(cbuf,state->x,0,state->natoms);
1707 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1709 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1710 bUpdateDoLR,fr->f_twin,fcd,
1711 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1712 wallcycle_stop(wcycle,ewcUPDATE);
1714 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1715 /* are the small terms in the shake_vir here due
1716 * to numerical errors, or are they important
1717 * physically? I'm thinking they are just errors, but not completely sure.
1718 * For now, will call without actually constraining, constr=NULL*/
1719 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1720 state,fr->bMolPBC,graph,f,
1721 &top->idef,tmp_vir,force_vir,
1722 cr,nrnb,wcycle,upd,NULL,
1723 bInitStep,FALSE,bCalcVir,
1726 if (!bOK && !bFFscan)
1728 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1731 if (fr->bSepDVDL && fplog && do_log)
1733 fprintf(fplog,sepdvdlformat,"Constraint dV/dl",0.0,dvdl);
1735 enerd->term[F_DVDL_BONDED] += dvdl;
1739 /* Need to unshift here */
1740 unshift_self(graph,state->box,state->x);
1743 GMX_BARRIER(cr->mpi_comm_mygroup);
1744 GMX_MPE_LOG(ev_update_finish);
1748 wallcycle_start(wcycle,ewcVSITECONSTR);
1751 shift_self(graph,state->box,state->x);
1753 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
1754 top->idef.iparams,top->idef.il,
1755 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
1759 unshift_self(graph,state->box,state->x);
1761 wallcycle_stop(wcycle,ewcVSITECONSTR);
1764 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1765 /* With Leap-Frog we can skip compute_globals at
1766 * non-communication steps, but we need to calculate
1767 * the kinetic energy one step before communication.
1769 if (bGStat || do_per_step(step+1,nstglobalcomm) ||
1772 if (ir->nstlist == -1 && bFirstIterate)
1774 gs.sig[eglsNABNSB] = nlh.nabnsb;
1776 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1777 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1779 bFirstIterate ? &gs : NULL,
1780 (step_rel % gs.nstms == 0) &&
1781 (multisim_nsteps<0 || (step_rel<multisim_nsteps)),
1783 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1785 | (!EI_VV(ir->eI) ? CGLO_ENERGY : 0)
1786 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1787 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1788 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1789 | (bIterations && iterate.bIterate ? CGLO_ITERATE : 0)
1790 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1793 if (ir->nstlist == -1 && bFirstIterate)
1795 nlh.nabnsb = gs.set[eglsNABNSB];
1796 gs.set[eglsNABNSB] = 0;
1799 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1800 /* ############# END CALC EKIN AND PRESSURE ################# */
1802 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1803 the virial that should probably be addressed eventually. state->veta has better properies,
1804 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1805 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1808 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1809 trace(shake_vir),&tracevir))
1813 bFirstIterate = FALSE;
1816 /* only add constraint dvdl after constraints */
1817 enerd->term[F_DVDL_BONDED] += dvdl;
1820 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1821 sum_dhdl(enerd,state->lambda,ir->fepvals);
1823 update_box(fplog,step,ir,mdatoms,state,graph,f,
1824 ir->nstlist==-1 ? &nlh.scale_tot : NULL,pcoupl_mu,nrnb,wcycle,upd,bInitStep,FALSE);
1826 /* ################# END UPDATE STEP 2 ################# */
1827 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1829 /* The coordinates (x) were unshifted in update */
1830 if (bFFscan && (shellfc==NULL || bConverged))
1832 if (print_forcefield(fplog,enerd->term,mdatoms->homenr,
1834 &(top_global->mols),mdatoms->massT,pres))
1838 fprintf(stderr,"\n");
1844 /* We will not sum ekinh_old,
1845 * so signal that we still have to do it.
1847 bSumEkinhOld = TRUE;
1852 /* Only do GCT when the relaxation of shells (minimization) has converged,
1853 * otherwise we might be coupling to bogus energies.
1854 * In parallel we must always do this, because the other sims might
1858 /* Since this is called with the new coordinates state->x, I assume
1859 * we want the new box state->box too. / EL 20040121
1861 do_coupling(fplog,oenv,nfile,fnm,tcr,t,step,enerd->term,fr,
1863 mdatoms,&(top->idef),mu_aver,
1864 top_global->mols.nr,cr,
1865 state->box,total_vir,pres,
1866 mu_tot,state->x,f,bConverged);
1870 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1872 /* use the directly determined last velocity, not actually the averaged half steps */
1873 if (bTrotter && ir->eI==eiVV)
1875 enerd->term[F_EKIN] = last_ekin;
1877 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1881 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1885 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir,state,&MassQ);
1887 /* Check for excessively large energies */
1891 real etot_max = 1e200;
1893 real etot_max = 1e30;
1895 if (fabs(enerd->term[F_ETOT]) > etot_max)
1897 fprintf(stderr,"Energy too large (%g), giving up\n",
1898 enerd->term[F_ETOT]);
1901 /* ######### END PREPARING EDR OUTPUT ########### */
1903 /* Time for performance */
1904 if (((step % stepout) == 0) || bLastStep)
1906 runtime_upd_proc(runtime);
1912 gmx_bool do_dr,do_or;
1914 if (fplog && do_log && bDoExpanded)
1916 /* only needed if doing expanded ensemble */
1917 PrintFreeEnergyInfoToFile(fplog,ir->fepvals,ir->expandedvals,ir->bSimTemp?ir->simtempvals:NULL,
1918 &df_history,state->fep_state,ir->nstlog,step);
1920 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1924 upd_mdebin(mdebin,bDoDHDL, TRUE,
1925 t,mdatoms->tmass,enerd,state,
1926 ir->fepvals,ir->expandedvals,lastbox,
1927 shake_vir,force_vir,total_vir,pres,
1928 ekind,mu_tot,constr);
1932 upd_mdebin_step(mdebin);
1935 do_dr = do_per_step(step,ir->nstdisreout);
1936 do_or = do_per_step(step,ir->nstorireout);
1938 print_ebin(outf->fp_ene,do_ene,do_dr,do_or,do_log?fplog:NULL,
1940 eprNORMAL,bCompact,mdebin,fcd,groups,&(ir->opts));
1942 if (ir->ePull != epullNO)
1944 pull_print_output(ir->pull,step,t);
1947 if (do_per_step(step,ir->nstlog))
1949 if(fflush(fplog) != 0)
1951 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of disk space?");
1957 /* Have to do this part after outputting the logfile and the edr file */
1958 state->fep_state = lamnew;
1959 for (i=0;i<efptNR;i++)
1961 state->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
1964 /* Remaining runtime */
1965 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1969 fprintf(stderr,"\n");
1971 print_time(stderr,runtime,step,ir,cr);
1974 /* Replica exchange */
1976 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1977 do_per_step(step,repl_ex_nst))
1979 bExchanged = replica_exchange(fplog,cr,repl_ex,
1983 if (bExchanged && DOMAINDECOMP(cr))
1985 dd_partition_system(fplog,step,cr,TRUE,1,
1986 state_global,top_global,ir,
1987 state,&f,mdatoms,top,fr,
1988 vsite,shellfc,constr,
1995 bStartingFromCpt = FALSE;
1997 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1998 /* With all integrators, except VV, we need to retain the pressure
1999 * at the current step for coupling at the next step.
2001 if ((state->flags & (1<<estPRES_PREV)) &&
2003 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
2005 /* Store the pressure in t_state for pressure coupling
2006 * at the next MD step.
2008 copy_mat(pres,state->pres_prev);
2011 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2013 if ( (membed!=NULL) && (!bLastStep) )
2015 rescale_membed(step_rel,membed,state_global->x);
2022 /* read next frame from input trajectory */
2023 bNotLastFrame = read_next_frame(oenv,status,&rerun_fr);
2028 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
2032 if (!bRerunMD || !rerun_fr.bStep)
2034 /* increase the MD step number */
2039 cycles = wallcycle_stop(wcycle,ewcSTEP);
2040 if (DOMAINDECOMP(cr) && wcycle)
2042 dd_cycles_add(cr->dd,cycles,ddCyclStep);
2045 if (bPMETuneRunning || bPMETuneTry)
2047 /* PME grid + cut-off optimization with GPUs or PME nodes */
2049 /* Count the total cycles over the last steps */
2050 cycles_pmes += cycles;
2052 /* We can only switch cut-off at NS steps */
2053 if (step % ir->nstlist == 0)
2055 /* PME grid + cut-off optimization with GPUs or PME nodes */
2058 if (DDMASTER(cr->dd))
2060 /* PME node load is too high, start tuning */
2061 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2063 dd_bcast(cr->dd,sizeof(gmx_bool),&bPMETuneRunning);
2065 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2067 bPMETuneTry = FALSE;
2070 if (bPMETuneRunning)
2072 /* init_step might not be a multiple of nstlist,
2073 * but the first cycle is always skipped anyhow.
2076 pme_load_balance(pme_loadbal,cr,
2077 (bVerbose && MASTER(cr)) ? stderr : NULL,
2079 ir,state,cycles_pmes,
2080 fr->ic,fr->nbv,&fr->pmedata,
2083 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2084 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2085 fr->rlist = fr->ic->rlist;
2086 fr->rlistlong = fr->ic->rlistlong;
2087 fr->rcoulomb = fr->ic->rcoulomb;
2088 fr->rvdw = fr->ic->rvdw;
2094 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2095 gs.set[eglsRESETCOUNTERS] != 0)
2097 /* Reset all the counters related to performance over the run */
2098 reset_all_counters(fplog,cr,step,&step_rel,ir,wcycle,nrnb,runtime,
2099 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2100 wcycle_set_reset_counters(wcycle,-1);
2101 /* Correct max_hours for the elapsed time */
2102 max_hours -= run_time/(60.0*60.0);
2103 bResetCountersHalfMaxH = FALSE;
2104 gs.set[eglsRESETCOUNTERS] = 0;
2108 /* End of main MD loop */
2112 runtime_end(runtime);
2114 if (bRerunMD && MASTER(cr))
2119 if (!(cr->duty & DUTY_PME))
2121 /* Tell the PME only node to finish */
2122 gmx_pme_send_finish(cr);
2127 if (ir->nstcalcenergy > 0 && !bRerunMD)
2129 print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
2130 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
2138 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2140 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)));
2141 fprintf(fplog,"Average number of atoms that crossed the half buffer length: %.1f\n\n",nlh.ab/nlh.nns);
2144 if (pme_loadbal != NULL)
2146 pme_loadbal_done(pme_loadbal,fplog);
2149 if (shellfc && fplog)
2151 fprintf(fplog,"Fraction of iterations that converged: %.2f %%\n",
2152 (nconverged*100.0)/step_rel);
2153 fprintf(fplog,"Average number of force evaluations per MD step: %.2f\n\n",
2157 if (repl_ex_nst > 0 && MASTER(cr))
2159 print_replica_exchange_statistics(fplog,repl_ex);
2162 runtime->nsteps_done = step_rel;