1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
56 #include "md_support.h"
57 #include "md_logging.h"
73 #include "domdec_network.h"
79 #include "compute_io.h"
81 #include "checkpoint.h"
82 #include "mtop_util.h"
83 #include "sighandler.h"
86 #include "pme_loadbal.h"
89 #include "types/nlistheuristics.h"
90 #include "types/iteratedconstraints.h"
91 #include "nbnxn_cuda_data_mgmt.h"
101 #include "corewrap.h"
104 static void reset_all_counters(FILE *fplog,t_commrec *cr,
105 gmx_large_int_t step,
106 gmx_large_int_t *step_rel,t_inputrec *ir,
107 gmx_wallcycle_t wcycle,t_nrnb *nrnb,
108 gmx_runtime_t *runtime,
109 nbnxn_cuda_ptr_t cu_nbv)
111 char sbuf[STEPSTRSIZE];
113 /* Reset all the counters related to performance over the run */
114 md_print_warn(cr,fplog,"step %s: resetting all time and cycle counters\n",
115 gmx_step_str(step,sbuf));
119 nbnxn_cuda_reset_timings(cu_nbv);
122 wallcycle_stop(wcycle,ewcRUN);
123 wallcycle_reset_all(wcycle);
124 if (DOMAINDECOMP(cr))
126 reset_dd_statistics_counters(cr->dd);
129 ir->init_step += *step_rel;
130 ir->nsteps -= *step_rel;
132 wallcycle_start(wcycle,ewcRUN);
133 runtime_start(runtime);
134 print_date_and_time(fplog,cr->nodeid,"Restarted time",runtime);
137 double do_md(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
138 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
140 gmx_vsite_t *vsite,gmx_constr_t constr,
141 int stepout,t_inputrec *ir,
142 gmx_mtop_t *top_global,
144 t_state *state_global,
146 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
147 gmx_edsam_t ed,t_forcerec *fr,
148 int repl_ex_nst,int repl_ex_nex,int repl_ex_seed,gmx_membed_t membed,
149 real cpt_period,real max_hours,
150 const char *deviceOptions,
152 gmx_runtime_t *runtime)
155 gmx_large_int_t step,step_rel;
157 double t,t0,lam0[efptNR];
158 gmx_bool bGStatEveryStep,bGStat,bCalcVir,bCalcEner;
159 gmx_bool bNS,bNStList,bSimAnn,bStopCM,bRerunMD,bNotLastFrame=FALSE,
160 bFirstStep,bStateFromCP,bStateFromTPX,bInitStep,bLastStep,
161 bBornRadii,bStartingFromCpt;
162 gmx_bool bDoDHDL=FALSE,bDoFEP=FALSE,bDoExpanded=FALSE;
163 gmx_bool do_ene,do_log,do_verbose,bRerunWarnNoV=TRUE,
164 bForceUpdate=FALSE,bCPT;
166 gmx_bool bMasterState;
167 int force_flags,cglo_flags;
168 tensor force_vir,shake_vir,total_vir,tmp_vir,pres;
173 t_state *bufstate=NULL;
174 matrix *scale_tot,pcoupl_mu,M,ebox;
177 gmx_repl_ex_t repl_ex=NULL;
180 t_mdebin *mdebin=NULL;
181 df_history_t df_history;
186 gmx_enerdata_t *enerd;
188 gmx_global_stat_t gstat;
189 gmx_update_t upd=NULL;
192 gmx_rng_t mcrng=NULL;
194 gmx_groups_t *groups;
195 gmx_ekindata_t *ekind, *ekind_save;
196 gmx_shellfc_t shellfc;
197 int count,nconverged=0;
200 gmx_bool bIonize=FALSE;
201 gmx_bool bTCR=FALSE,bConverged=TRUE,bOK,bSumEkinhOld,bExchanged;
203 gmx_bool bResetCountersHalfMaxH=FALSE;
204 gmx_bool bVV,bIterations,bFirstIterate,bTemp,bPres,bTrotter;
205 gmx_bool bUpdateDoLR;
208 atom_id *grpindex=NULL;
210 t_coupl_rec *tcr=NULL;
211 rvec *xcopy=NULL,*vcopy=NULL,*cbuf=NULL;
212 matrix boxcopy={{0}},lastbox;
214 real fom,oldfom,veta_save,pcurr,scalevir,tracevir;
222 real saved_conserved_quantity = 0;
227 char sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
228 int handled_stop_condition=gmx_stop_cond_none; /* compare to get_stop_condition*/
229 gmx_iterate_t iterate;
230 gmx_large_int_t multisim_nsteps=-1; /* number of steps to do before first multisim
231 simulation stops. If equal to zero, don't
232 communicate any more between multisims.*/
233 /* PME load balancing data for GPU kernels */
234 pme_load_balancing_t pme_loadbal=NULL;
236 gmx_bool bPMETuneTry=FALSE,bPMETuneRunning=FALSE;
239 /* Temporary addition for FAHCORE checkpointing */
243 /* Check for special mdrun options */
244 bRerunMD = (Flags & MD_RERUN);
245 bIonize = (Flags & MD_IONIZE);
246 bFFscan = (Flags & MD_FFSCAN);
247 bAppend = (Flags & MD_APPENDFILES);
248 if (Flags & MD_RESETCOUNTERSHALFWAY)
252 /* Signal to reset the counters half the simulation steps. */
253 wcycle_set_reset_counters(wcycle,ir->nsteps/2);
255 /* Signal to reset the counters halfway the simulation time. */
256 bResetCountersHalfMaxH = (max_hours > 0);
259 /* md-vv uses averaged full step velocities for T-control
260 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
261 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
263 if (bVV) /* to store the initial velocities while computing virial */
265 snew(cbuf,top_global->natoms);
267 /* all the iteratative cases - only if there are constraints */
268 bIterations = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
269 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
273 /* Since we don't know if the frames read are related in any way,
274 * rebuild the neighborlist at every step.
277 ir->nstcalcenergy = 1;
281 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
283 nstglobalcomm = check_nstglobalcomm(fplog,cr,nstglobalcomm,ir);
284 bGStatEveryStep = (nstglobalcomm == 1);
286 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
289 "To reduce the energy communication with nstlist = -1\n"
290 "the neighbor list validity should not be checked at every step,\n"
291 "this means that exact integration is not guaranteed.\n"
292 "The neighbor list validity is checked after:\n"
293 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
294 "In most cases this will result in exact integration.\n"
295 "This reduces the energy communication by a factor of 2 to 3.\n"
296 "If you want less energy communication, set nstlist > 3.\n\n");
299 if (bRerunMD || bFFscan)
303 groups = &top_global->groups;
306 init_md(fplog,cr,ir,oenv,&t,&t0,state_global->lambda,
307 &(state_global->fep_state),lam0,
308 nrnb,top_global,&upd,
309 nfile,fnm,&outf,&mdebin,
310 force_vir,shake_vir,mu_tot,&bSimAnn,&vcm,state_global,Flags);
312 clear_mat(total_vir);
314 /* Energy terms and groups */
316 init_enerdata(top_global->groups.grps[egcENER].nr,ir->fepvals->n_lambda,
318 if (DOMAINDECOMP(cr))
324 snew(f,top_global->natoms);
327 /* lambda Monte carlo random number generator */
330 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
332 /* copy the state into df_history */
333 copy_df_history(&df_history,&state_global->dfhist);
335 /* Kinetic energy data */
337 init_ekindata(fplog,top_global,&(ir->opts),ekind);
338 /* needed for iteration of constraints */
340 init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
341 /* Copy the cos acceleration to the groups struct */
342 ekind->cosacc.cos_accel = ir->cos_accel;
344 gstat = global_stat_init(ir);
347 /* Check for polarizable models and flexible constraints */
348 shellfc = init_shell_flexcon(fplog,
349 top_global,n_flexible_constraints(constr),
350 (ir->bContinuation ||
351 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
352 NULL : state_global->x);
356 #ifdef GMX_THREAD_MPI
357 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
359 set_deform_reference_box(upd,
360 deform_init_init_step_tpx,
361 deform_init_box_tpx);
362 #ifdef GMX_THREAD_MPI
363 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
368 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
369 if ((io > 2000) && MASTER(cr))
371 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
375 if (DOMAINDECOMP(cr)) {
376 top = dd_init_local_top(top_global);
379 dd_init_local_state(cr->dd,state_global,state);
381 if (DDMASTER(cr->dd) && ir->nstfout) {
382 snew(f_global,state_global->natoms);
386 /* Initialize the particle decomposition and split the topology */
387 top = split_system(fplog,top_global,ir,cr);
389 pd_cg_range(cr,&fr->cg0,&fr->hcg);
390 pd_at_range(cr,&a0,&a1);
392 top = gmx_mtop_generate_local_top(top_global,ir);
395 a1 = top_global->natoms;
398 forcerec_set_excl_load(fr,top,cr);
400 state = partdec_init_local_state(cr,state_global);
403 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
406 set_vsite_top(vsite,top,mdatoms,cr);
409 if (ir->ePBC != epbcNONE && !fr->bMolPBC) {
410 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
414 make_local_shells(cr,mdatoms,shellfc);
417 init_bonded_thread_force_reduction(fr,&top->idef);
419 if (ir->pull && PAR(cr)) {
420 dd_make_local_pull_groups(NULL,ir->pull,mdatoms);
424 if (DOMAINDECOMP(cr))
426 /* Distribute the charge groups over the nodes from the master node */
427 dd_partition_system(fplog,ir->init_step,cr,TRUE,1,
428 state_global,top_global,ir,
429 state,&f,mdatoms,top,fr,
430 vsite,shellfc,constr,
435 update_mdatoms(mdatoms,state->lambda[efptMASS]);
437 if (opt2bSet("-cpi",nfile,fnm))
439 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr);
443 bStateFromCP = FALSE;
450 /* Update mdebin with energy history if appending to output files */
451 if ( Flags & MD_APPENDFILES )
453 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
457 /* We might have read an energy history from checkpoint,
458 * free the allocated memory and reset the counts.
460 done_energyhistory(&state_global->enerhist);
461 init_energyhistory(&state_global->enerhist);
464 /* Set the initial energy history in state by updating once */
465 update_energyhistory(&state_global->enerhist,mdebin);
468 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
470 /* Set the random state if we read a checkpoint file */
471 set_stochd_state(upd,state);
474 if (state->flags & (1<<estMC_RNG))
476 set_mc_state(mcrng,state);
479 /* Initialize constraints */
481 if (!DOMAINDECOMP(cr))
482 set_constraints(constr,top,ir,mdatoms,cr);
485 /* Check whether we have to GCT stuff */
486 bTCR = ftp2bSet(efGCT,nfile,fnm);
489 fprintf(stderr,"Will do General Coupling Theory!\n");
491 gnx = top_global->mols.nr;
493 for(i=0; (i<gnx); i++) {
500 /* We need to be sure replica exchange can only occur
501 * when the energies are current */
502 check_nst_param(fplog,cr,"nstcalcenergy",ir->nstcalcenergy,
503 "repl_ex_nst",&repl_ex_nst);
504 /* This check needs to happen before inter-simulation
505 * signals are initialized, too */
507 if (repl_ex_nst > 0 && MASTER(cr))
509 repl_ex = init_replica_exchange(fplog,cr->ms,state_global,ir,
510 repl_ex_nst,repl_ex_nex,repl_ex_seed);
513 /* PME tuning is only supported with GPUs or PME nodes and not with rerun */
514 if ((Flags & MD_TUNEPME) &&
515 EEL_PME(fr->eeltype) &&
516 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
519 pme_loadbal_init(&pme_loadbal,ir,state->box,fr->ic,fr->pmedata);
521 if (cr->duty & DUTY_PME)
523 /* Start tuning right away, as we can't measure the load */
524 bPMETuneRunning = TRUE;
528 /* Separate PME nodes, we can measure the PP/PME load balance */
533 if (!ir->bContinuation && !bRerunMD)
535 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
537 /* Set the velocities of frozen particles to zero */
538 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
542 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
552 /* Constrain the initial coordinates and velocities */
553 do_constrain_first(fplog,constr,ir,mdatoms,state,f,
554 graph,cr,nrnb,fr,top,shake_vir);
558 /* Construct the virtual sites for the initial configuration */
559 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
560 top->idef.iparams,top->idef.il,
561 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
567 /* set free energy calculation frequency as the minimum of nstdhdl, nstexpanded, and nstrepl_ex_nst*/
568 nstfep = ir->fepvals->nstdhdl;
569 if (ir->bExpanded && (nstfep > ir->expandedvals->nstexpanded))
571 nstfep = ir->expandedvals->nstexpanded;
573 if (repl_ex_nst > 0 && nstfep > repl_ex_nst)
575 nstfep = repl_ex_nst;
578 /* I'm assuming we need global communication the first time! MRS */
579 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
580 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM:0)
581 | (bVV ? CGLO_PRESSURE:0)
582 | (bVV ? CGLO_CONSTRAINT:0)
583 | (bRerunMD ? CGLO_RERUNMD:0)
584 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN:0));
586 bSumEkinhOld = FALSE;
587 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
588 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
589 constr,NULL,FALSE,state->box,
590 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,cglo_flags);
591 if (ir->eI == eiVVAK) {
592 /* a second call to get the half step temperature initialized as well */
593 /* we do the same call as above, but turn the pressure off -- internally to
594 compute_globals, this is recognized as a velocity verlet half-step
595 kinetic energy calculation. This minimized excess variables, but
596 perhaps loses some logic?*/
598 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
599 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
600 constr,NULL,FALSE,state->box,
601 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
602 cglo_flags &~ (CGLO_STOPCM | CGLO_PRESSURE));
605 /* Calculate the initial half step temperature, and save the ekinh_old */
606 if (!(Flags & MD_STARTFROMCPT))
608 for(i=0; (i<ir->opts.ngtc); i++)
610 copy_mat(ekind->tcstat[i].ekinh,ekind->tcstat[i].ekinh_old);
615 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
616 and there is no previous step */
619 /* if using an iterative algorithm, we need to create a working directory for the state. */
622 bufstate = init_bufstate(state);
626 snew(xcopy,state->natoms);
627 snew(vcopy,state->natoms);
628 copy_rvecn(state->x,xcopy,0,state->natoms);
629 copy_rvecn(state->v,vcopy,0,state->natoms);
630 copy_mat(state->box,boxcopy);
633 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
634 temperature control */
635 trotter_seq = init_npt_vars(ir,state,&MassQ,bTrotter);
639 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
642 "RMS relative constraint deviation after constraining: %.2e\n",
643 constr_rmsd(constr,FALSE));
645 if (EI_STATE_VELOCITY(ir->eI))
647 fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
651 fprintf(stderr,"starting md rerun '%s', reading coordinates from"
652 " input trajectory '%s'\n\n",
653 *(top_global->name),opt2fn("-rerun",nfile,fnm));
656 fprintf(stderr,"Calculated time to finish depends on nsteps from "
657 "run input file,\nwhich may not correspond to the time "
658 "needed to process input trajectory.\n\n");
664 fprintf(stderr,"starting mdrun '%s'\n",
665 *(top_global->name));
668 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
672 sprintf(tbuf,"%s","infinite");
674 if (ir->init_step > 0)
676 fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
677 gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
678 gmx_step_str(ir->init_step,sbuf2),
679 ir->init_step*ir->delta_t);
683 fprintf(stderr,"%s steps, %s ps.\n",
684 gmx_step_str(ir->nsteps,sbuf),tbuf);
690 /* Set and write start time */
691 runtime_start(runtime);
692 print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
693 wallcycle_start(wcycle,ewcRUN);
699 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
701 chkpt_ret=fcCheckPointParallel( cr->nodeid,
703 if ( chkpt_ret == 0 )
704 gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", 0 );
708 /***********************************************************
712 ************************************************************/
714 /* if rerunMD then read coordinates and velocities from input trajectory */
717 if (getenv("GMX_FORCE_UPDATE"))
725 bNotLastFrame = read_first_frame(oenv,&status,
726 opt2fn("-rerun",nfile,fnm),
727 &rerun_fr,TRX_NEED_X | TRX_READ_V);
728 if (rerun_fr.natoms != top_global->natoms)
731 "Number of atoms in trajectory (%d) does not match the "
732 "run input file (%d)\n",
733 rerun_fr.natoms,top_global->natoms);
735 if (ir->ePBC != epbcNONE)
739 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);
741 if (max_cutoff2(ir->ePBC,rerun_fr.box) < sqr(fr->rlistlong))
743 gmx_fatal(FARGS,"Rerun trajectory frame step %d time %f has too small box dimensions",rerun_fr.step,rerun_fr.time);
750 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
753 if (ir->ePBC != epbcNONE)
755 /* Set the shift vectors.
756 * Necessary here when have a static box different from the tpr box.
758 calc_shifts(rerun_fr.box,fr->shift_vec);
762 /* loop over MD steps or if rerunMD to end of input trajectory */
764 /* Skip the first Nose-Hoover integration when we get the state from tpx */
765 bStateFromTPX = !bStateFromCP;
766 bInitStep = bFirstStep && (bStateFromTPX || bVV);
767 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
769 bSumEkinhOld = FALSE;
772 init_global_signals(&gs,cr,ir,repl_ex_nst);
774 step = ir->init_step;
777 if (ir->nstlist == -1)
779 init_nlistheuristics(&nlh,bGStatEveryStep,step);
782 if (MULTISIM(cr) && (repl_ex_nst <=0 ))
784 /* check how many steps are left in other sims */
785 multisim_nsteps=get_multisim_nsteps(cr, ir->nsteps);
789 /* and stop now if we should */
790 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
791 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
792 while (!bLastStep || (bRerunMD && bNotLastFrame)) {
794 wallcycle_start(wcycle,ewcSTEP);
797 if (rerun_fr.bStep) {
798 step = rerun_fr.step;
799 step_rel = step - ir->init_step;
801 if (rerun_fr.bTime) {
811 bLastStep = (step_rel == ir->nsteps);
812 t = t0 + step*ir->delta_t;
815 if (ir->efep != efepNO || ir->bSimTemp)
817 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
818 requiring different logic. */
820 set_current_lambdas(step,ir->fepvals,bRerunMD,&rerun_fr,state_global,state,lam0);
821 bDoDHDL = do_per_step(step,ir->fepvals->nstdhdl);
822 bDoFEP = (do_per_step(step,nstfep) && (ir->efep != efepNO));
823 bDoExpanded = (do_per_step(step,ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
828 update_annealing_target_temp(&(ir->opts),t);
833 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
835 for(i=0; i<state_global->natoms; i++)
837 copy_rvec(rerun_fr.x[i],state_global->x[i]);
841 for(i=0; i<state_global->natoms; i++)
843 copy_rvec(rerun_fr.v[i],state_global->v[i]);
848 for(i=0; i<state_global->natoms; i++)
850 clear_rvec(state_global->v[i]);
854 fprintf(stderr,"\nWARNING: Some frames do not contain velocities.\n"
855 " Ekin, temperature and pressure are incorrect,\n"
856 " the virial will be incorrect when constraints are present.\n"
858 bRerunWarnNoV = FALSE;
862 copy_mat(rerun_fr.box,state_global->box);
863 copy_mat(state_global->box,state->box);
865 if (vsite && (Flags & MD_RERUN_VSITE))
867 if (DOMAINDECOMP(cr))
869 gmx_fatal(FARGS,"Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
873 /* Following is necessary because the graph may get out of sync
874 * with the coordinates if we only have every N'th coordinate set
876 mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
877 shift_self(graph,state->box,state->x);
879 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
880 top->idef.iparams,top->idef.il,
881 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
884 unshift_self(graph,state->box,state->x);
889 /* Stop Center of Mass motion */
890 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step,ir->nstcomm));
892 /* Copy back starting coordinates in case we're doing a forcefield scan */
895 for(ii=0; (ii<state->natoms); ii++)
897 copy_rvec(xcopy[ii],state->x[ii]);
898 copy_rvec(vcopy[ii],state->v[ii]);
900 copy_mat(boxcopy,state->box);
905 /* for rerun MD always do Neighbour Searching */
906 bNS = (bFirstStep || ir->nstlist != 0);
911 /* Determine whether or not to do Neighbour Searching and LR */
912 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
914 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
915 (ir->nstlist == -1 && nlh.nabnsb > 0));
917 if (bNS && ir->nstlist == -1)
919 set_nlistheuristics(&nlh,bFirstStep || bExchanged || bDoFEP, step);
923 /* check whether we should stop because another simulation has
927 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
928 (multisim_nsteps != ir->nsteps) )
935 "Stopping simulation %d because another one has finished\n",
939 gs.sig[eglsCHKPT] = 1;
944 /* < 0 means stop at next step, > 0 means stop at next NS step */
945 if ( (gs.set[eglsSTOPCOND] < 0 ) ||
946 ( (gs.set[eglsSTOPCOND] > 0 ) && ( bNS || ir->nstlist==0)) )
951 /* Determine whether or not to update the Born radii if doing GB */
952 bBornRadii=bFirstStep;
953 if (ir->implicit_solvent && (step % ir->nstgbradii==0))
958 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
959 do_verbose = bVerbose &&
960 (step % stepout == 0 || bFirstStep || bLastStep);
962 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
970 bMasterState = FALSE;
971 /* Correct the new box if it is too skewed */
972 if (DYNAMIC_BOX(*ir))
974 if (correct_box(fplog,step,state->box,graph))
979 if (DOMAINDECOMP(cr) && bMasterState)
981 dd_collect_state(cr->dd,state,state_global);
985 if (DOMAINDECOMP(cr))
987 /* Repartition the domain decomposition */
988 wallcycle_start(wcycle,ewcDOMDEC);
989 dd_partition_system(fplog,step,cr,
990 bMasterState,nstglobalcomm,
991 state_global,top_global,ir,
992 state,&f,mdatoms,top,fr,
993 vsite,shellfc,constr,
995 do_verbose && !bPMETuneRunning);
996 wallcycle_stop(wcycle,ewcDOMDEC);
997 /* If using an iterative integrator, reallocate space to match the decomposition */
1001 if (MASTER(cr) && do_log && !bFFscan)
1003 print_ebin_header(fplog,step,t,state->lambda[efptFEP]); /* can we improve the information printed here? */
1006 if (ir->efep != efepNO)
1008 update_mdatoms(mdatoms,state->lambda[efptMASS]);
1011 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1014 /* We need the kinetic energy at minus the half step for determining
1015 * the full step kinetic energy and possibly for T-coupling.*/
1016 /* This may not be quite working correctly yet . . . . */
1017 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1018 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1019 constr,NULL,FALSE,state->box,
1020 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1021 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1023 clear_mat(force_vir);
1025 /* Ionize the atoms if necessary */
1028 ionize(fplog,oenv,mdatoms,top_global,t,ir,state->x,state->v,
1029 mdatoms->start,mdatoms->start+mdatoms->homenr,state->box,cr);
1032 /* Update force field in ffscan program */
1035 if (update_forcefield(fplog,
1037 mdatoms->nr,state->x,state->box))
1045 /* We write a checkpoint at this MD step when:
1046 * either at an NS step when we signalled through gs,
1047 * or at the last step (but not when we do not want confout),
1048 * but never at the first step or with rerun.
1050 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1051 (bLastStep && (Flags & MD_CONFOUT))) &&
1052 step > ir->init_step && !bRerunMD);
1055 gs.set[eglsCHKPT] = 0;
1058 /* Determine the energy and pressure:
1059 * at nstcalcenergy steps and at energy output steps (set below).
1061 if (EI_VV(ir->eI) && (!bInitStep))
1063 /* for vv, the first half actually corresponds to the last step */
1064 bCalcEner = do_per_step(step-1,ir->nstcalcenergy);
1068 bCalcEner = do_per_step(step,ir->nstcalcenergy);
1070 bCalcVir = bCalcEner ||
1071 (ir->epc != epcNO && do_per_step(step,ir->nstpcouple));
1073 /* Do we need global communication ? */
1074 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1075 do_per_step(step,nstglobalcomm) ||
1076 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1078 do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
1080 if (do_ene || do_log)
1087 /* these CGLO_ options remain the same throughout the iteration */
1088 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1089 (bGStat ? CGLO_GSTAT : 0)
1092 force_flags = (GMX_FORCE_STATECHANGED |
1093 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1094 GMX_FORCE_ALLFORCES |
1096 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1097 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1098 (bDoFEP ? GMX_FORCE_DHDL : 0)
1103 if(do_per_step(step,ir->nstcalclr))
1105 force_flags |= GMX_FORCE_DO_LR;
1111 /* Now is the time to relax the shells */
1112 count=relax_shell_flexcon(fplog,cr,bVerbose,bFFscan ? step+1 : step,
1114 bStopCM,top,top_global,
1116 state,f,force_vir,mdatoms,
1117 nrnb,wcycle,graph,groups,
1118 shellfc,fr,bBornRadii,t,mu_tot,
1119 state->natoms,&bConverged,vsite,
1130 /* The coordinates (x) are shifted (to get whole molecules)
1132 * This is parallellized as well, and does communication too.
1133 * Check comments in sim_util.c
1135 do_force(fplog,cr,ir,step,nrnb,wcycle,top,top_global,groups,
1136 state->box,state->x,&state->hist,
1137 f,force_vir,mdatoms,enerd,fcd,
1138 state->lambda,graph,
1139 fr,vsite,mu_tot,t,outf->fp_field,ed,bBornRadii,
1140 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1145 mu_aver = calc_mu_aver(cr,state->x,mdatoms->chargeA,
1146 mu_tot,&top_global->mols,mdatoms,gnx,grpindex);
1149 if (bTCR && bFirstStep)
1151 tcr=init_coupling(fplog,nfile,fnm,cr,fr,mdatoms,&(top->idef));
1152 fprintf(fplog,"Done init_coupling\n");
1156 if (bVV && !bStartingFromCpt && !bRerunMD)
1157 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1159 if (ir->eI==eiVV && bInitStep)
1161 /* if using velocity verlet with full time step Ekin,
1162 * take the first half step only to compute the
1163 * virial for the first step. From there,
1164 * revert back to the initial coordinates
1165 * so that the input is actually the initial step.
1167 copy_rvecn(state->v,cbuf,0,state->natoms); /* should make this better for parallelizing? */
1169 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1170 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ1);
1173 /* If we are using twin-range interactions where the long-range component
1174 * is only evaluated every nstcalclr>1 steps, we should do a special update
1175 * step to combine the long-range forces on these steps.
1176 * For nstcalclr=1 this is not done, since the forces would have been added
1177 * directly to the short-range forces already.
1179 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1181 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,
1182 f,bUpdateDoLR,fr->f_twin,fcd,
1183 ekind,M,wcycle,upd,bInitStep,etrtVELOCITY1,
1184 cr,nrnb,constr,&top->idef);
1188 gmx_iterate_init(&iterate,bIterations && !bInitStep);
1190 /* for iterations, we save these vectors, as we will be self-consistently iterating
1193 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1195 /* save the state */
1196 if (bIterations && iterate.bIterate) {
1197 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1200 bFirstIterate = TRUE;
1201 while (bFirstIterate || (bIterations && iterate.bIterate))
1203 if (bIterations && iterate.bIterate)
1205 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1206 if (bFirstIterate && bTrotter)
1208 /* The first time through, we need a decent first estimate
1209 of veta(t+dt) to compute the constraints. Do
1210 this by computing the box volume part of the
1211 trotter integration at this time. Nothing else
1212 should be changed by this routine here. If
1213 !(first time), we start with the previous value
1216 veta_save = state->veta;
1217 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ0);
1218 vetanew = state->veta;
1219 state->veta = veta_save;
1224 if ( !bRerunMD || rerun_fr.bV || bForceUpdate) { /* Why is rerun_fr.bV here? Unclear. */
1227 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1228 state,fr->bMolPBC,graph,f,
1229 &top->idef,shake_vir,NULL,
1230 cr,nrnb,wcycle,upd,constr,
1231 bInitStep,TRUE,bCalcVir,vetanew);
1233 if (!bOK && !bFFscan)
1235 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1240 { /* Need to unshift here if a do_force has been
1241 called in the previous step */
1242 unshift_self(graph,state->box,state->x);
1246 /* if VV, compute the pressure and constraints */
1247 /* For VV2, we strictly only need this if using pressure
1248 * control, but we really would like to have accurate pressures
1250 * Think about ways around this in the future?
1251 * For now, keep this choice in comments.
1253 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1254 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1256 bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK));
1257 if (bCalcEner && ir->eI==eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1259 bSumEkinhOld = TRUE;
1261 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1262 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1263 constr,NULL,FALSE,state->box,
1264 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1267 | (bTemp ? CGLO_TEMPERATURE:0)
1268 | (bPres ? CGLO_PRESSURE : 0)
1269 | (bPres ? CGLO_CONSTRAINT : 0)
1270 | ((bIterations && iterate.bIterate) ? CGLO_ITERATE : 0)
1271 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1274 /* explanation of above:
1275 a) We compute Ekin at the full time step
1276 if 1) we are using the AveVel Ekin, and it's not the
1277 initial step, or 2) if we are using AveEkin, but need the full
1278 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1279 b) If we are using EkinAveEkin for the kinetic energy for the temperture control, we still feed in
1280 EkinAveVel because it's needed for the pressure */
1282 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1287 m_add(force_vir,shake_vir,total_vir); /* we need the un-dispersion corrected total vir here */
1288 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ2);
1295 /* We need the kinetic energy at minus the half step for determining
1296 * the full step kinetic energy and possibly for T-coupling.*/
1297 /* This may not be quite working correctly yet . . . . */
1298 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1299 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1300 constr,NULL,FALSE,state->box,
1301 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1302 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1308 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1309 state->veta,&vetanew))
1313 bFirstIterate = FALSE;
1316 if (bTrotter && !bInitStep) {
1317 enerd->term[F_DVDL_BONDED] += dvdl; /* only add after iterations */
1318 copy_mat(shake_vir,state->svir_prev);
1319 copy_mat(force_vir,state->fvir_prev);
1320 if (IR_NVT_TROTTER(ir) && ir->eI==eiVV) {
1321 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1322 enerd->term[F_TEMP] = sum_ekin(&(ir->opts),ekind,NULL,(ir->eI==eiVV),FALSE,FALSE);
1323 enerd->term[F_EKIN] = trace(ekind->ekin);
1326 /* if it's the initial step, we performed this first step just to get the constraint virial */
1327 if (bInitStep && ir->eI==eiVV) {
1328 copy_rvecn(cbuf,state->v,0,state->natoms);
1331 if (fr->bSepDVDL && fplog && do_log)
1333 fprintf(fplog,sepdvdlformat,"Constraint",0.0,dvdl);
1335 enerd->term[F_DVDL_BONDED] += dvdl;
1338 /* MRS -- now done iterating -- compute the conserved quantity */
1340 saved_conserved_quantity = compute_conserved_from_auxiliary(ir,state,&MassQ);
1343 last_ekin = enerd->term[F_EKIN];
1345 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1347 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1349 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1352 sum_dhdl(enerd,state->lambda,ir->fepvals);
1356 /* ######## END FIRST UPDATE STEP ############## */
1357 /* ######## If doing VV, we now have v(dt) ###### */
1359 /* perform extended ensemble sampling in lambda - we don't
1360 actually move to the new state before outputting
1361 statistics, but if performing simulated tempering, we
1362 do update the velocities and the tau_t. */
1364 lamnew = ExpandedEnsembleDynamics(fplog,ir,enerd,state,&MassQ,&df_history,step,mcrng,state->v,mdatoms);
1366 /* ################## START TRAJECTORY OUTPUT ################# */
1368 /* Now we have the energies and forces corresponding to the
1369 * coordinates at time t. We must output all of this before
1371 * for RerunMD t is read from input trajectory
1374 if (do_per_step(step,ir->nstxout)) { mdof_flags |= MDOF_X; }
1375 if (do_per_step(step,ir->nstvout)) { mdof_flags |= MDOF_V; }
1376 if (do_per_step(step,ir->nstfout)) { mdof_flags |= MDOF_F; }
1377 if (do_per_step(step,ir->nstxtcout)) { mdof_flags |= MDOF_XTC; }
1378 if (bCPT) { mdof_flags |= MDOF_CPT; };
1380 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1383 /* Enforce writing positions and velocities at end of run */
1384 mdof_flags |= (MDOF_X | MDOF_V);
1389 fcReportProgress( ir->nsteps, step );
1391 /* sync bCPT and fc record-keeping */
1392 if (bCPT && MASTER(cr))
1393 fcRequestCheckPoint();
1396 if (mdof_flags != 0)
1398 wallcycle_start(wcycle,ewcTRAJ);
1401 if (state->flags & (1<<estLD_RNG))
1403 get_stochd_state(upd,state);
1405 if (state->flags & (1<<estMC_RNG))
1407 get_mc_state(mcrng,state);
1413 state_global->ekinstate.bUpToDate = FALSE;
1417 update_ekinstate(&state_global->ekinstate,ekind);
1418 state_global->ekinstate.bUpToDate = TRUE;
1420 update_energyhistory(&state_global->enerhist,mdebin);
1421 if (ir->efep!=efepNO || ir->bSimTemp)
1423 state_global->fep_state = state->fep_state; /* MRS: seems kludgy. The code should be
1424 structured so this isn't necessary.
1425 Note this reassignment is only necessary
1426 for single threads.*/
1427 copy_df_history(&state_global->dfhist,&df_history);
1431 write_traj(fplog,cr,outf,mdof_flags,top_global,
1432 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
1439 if (bLastStep && step_rel == ir->nsteps &&
1440 (Flags & MD_CONFOUT) && MASTER(cr) &&
1441 !bRerunMD && !bFFscan)
1443 /* x and v have been collected in write_traj,
1444 * because a checkpoint file will always be written
1447 fprintf(stderr,"\nWriting final coordinates.\n");
1450 /* Make molecules whole only for confout writing */
1451 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
1453 write_sto_conf_mtop(ftp2fn(efSTO,nfile,fnm),
1454 *top_global->name,top_global,
1455 state_global->x,state_global->v,
1456 ir->ePBC,state->box);
1459 wallcycle_stop(wcycle,ewcTRAJ);
1462 /* kludge -- virial is lost with restart for NPT control. Must restart */
1463 if (bStartingFromCpt && bVV)
1465 copy_mat(state->svir_prev,shake_vir);
1466 copy_mat(state->fvir_prev,force_vir);
1468 /* ################## END TRAJECTORY OUTPUT ################ */
1470 /* Determine the wallclock run time up till now */
1471 run_time = gmx_gettime() - (double)runtime->real;
1473 /* Check whether everything is still allright */
1474 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1475 #ifdef GMX_THREAD_MPI
1480 /* this is just make gs.sig compatible with the hack
1481 of sending signals around by MPI_Reduce with together with
1483 if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
1484 gs.sig[eglsSTOPCOND]=1;
1485 if ( gmx_get_stop_condition() == gmx_stop_cond_next )
1486 gs.sig[eglsSTOPCOND]=-1;
1487 /* < 0 means stop at next step, > 0 means stop at next NS step */
1491 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1492 gmx_get_signal_name(),
1493 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1497 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1498 gmx_get_signal_name(),
1499 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1501 handled_stop_condition=(int)gmx_get_stop_condition();
1503 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1504 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1505 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1507 /* Signal to terminate the run */
1508 gs.sig[eglsSTOPCOND] = 1;
1511 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1513 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1516 if (bResetCountersHalfMaxH && MASTER(cr) &&
1517 run_time > max_hours*60.0*60.0*0.495)
1519 gs.sig[eglsRESETCOUNTERS] = 1;
1522 if (ir->nstlist == -1 && !bRerunMD)
1524 /* When bGStatEveryStep=FALSE, global_stat is only called
1525 * when we check the atom displacements, not at NS steps.
1526 * This means that also the bonded interaction count check is not
1527 * performed immediately after NS. Therefore a few MD steps could
1528 * be performed with missing interactions.
1529 * But wrong energies are never written to file,
1530 * since energies are only written after global_stat
1533 if (step >= nlh.step_nscheck)
1535 nlh.nabnsb = natoms_beyond_ns_buffer(ir,fr,&top->cgs,
1536 nlh.scale_tot,state->x);
1540 /* This is not necessarily true,
1541 * but step_nscheck is determined quite conservatively.
1547 /* In parallel we only have to check for checkpointing in steps
1548 * where we do global communication,
1549 * otherwise the other nodes don't know.
1551 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1554 run_time >= nchkpt*cpt_period*60.0)) &&
1555 gs.set[eglsCHKPT] == 0)
1557 gs.sig[eglsCHKPT] = 1;
1560 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1565 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1567 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1569 gmx_bool bIfRandomize;
1570 bIfRandomize = update_randomize_velocities(ir,step,mdatoms,state,upd,&top->idef,constr);
1571 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1572 if (constr && bIfRandomize)
1574 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1575 state,fr->bMolPBC,graph,f,
1576 &top->idef,tmp_vir,NULL,
1577 cr,nrnb,wcycle,upd,constr,
1578 bInitStep,TRUE,bCalcVir,vetanew);
1585 gmx_iterate_init(&iterate,bIterations);
1588 /* for iterations, we save these vectors, as we will be redoing the calculations */
1589 if (bIterations && iterate.bIterate)
1591 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1593 bFirstIterate = TRUE;
1594 while (bFirstIterate || (bIterations && iterate.bIterate))
1596 /* We now restore these vectors to redo the calculation with improved extended variables */
1599 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1602 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1603 so scroll down for that logic */
1605 /* ######### START SECOND UPDATE STEP ################# */
1606 /* Box is changed in update() when we do pressure coupling,
1607 * but we should still use the old box for energy corrections and when
1608 * writing it to the energy file, so it matches the trajectory files for
1609 * the same timestep above. Make a copy in a separate array.
1611 copy_mat(state->box,lastbox);
1614 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1616 wallcycle_start(wcycle,ewcUPDATE);
1618 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1621 if (bIterations && iterate.bIterate)
1629 /* we use a new value of scalevir to converge the iterations faster */
1630 scalevir = tracevir/trace(shake_vir);
1632 msmul(shake_vir,scalevir,shake_vir);
1633 m_add(force_vir,shake_vir,total_vir);
1634 clear_mat(shake_vir);
1636 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ3);
1637 /* We can only do Berendsen coupling after we have summed
1638 * the kinetic energy or virial. Since the happens
1639 * in global_state after update, we should only do it at
1640 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1645 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1646 update_pcouple(fplog,step,ir,state,pcoupl_mu,M,wcycle,
1652 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1654 /* velocity half-step update */
1655 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1656 bUpdateDoLR,fr->f_twin,fcd,
1657 ekind,M,wcycle,upd,FALSE,etrtVELOCITY2,
1658 cr,nrnb,constr,&top->idef);
1661 /* Above, initialize just copies ekinh into ekin,
1662 * it doesn't copy position (for VV),
1663 * and entire integrator for MD.
1668 copy_rvecn(state->x,cbuf,0,state->natoms);
1670 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1672 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1673 bUpdateDoLR,fr->f_twin,fcd,
1674 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1675 wallcycle_stop(wcycle,ewcUPDATE);
1677 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,state,
1678 fr->bMolPBC,graph,f,
1679 &top->idef,shake_vir,force_vir,
1680 cr,nrnb,wcycle,upd,constr,
1681 bInitStep,FALSE,bCalcVir,state->veta);
1685 /* erase F_EKIN and F_TEMP here? */
1686 /* just compute the kinetic energy at the half step to perform a trotter step */
1687 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1688 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1689 constr,NULL,FALSE,lastbox,
1690 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1691 cglo_flags | CGLO_TEMPERATURE
1693 wallcycle_start(wcycle,ewcUPDATE);
1694 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ4);
1695 /* now we know the scaling, we can compute the positions again again */
1696 copy_rvecn(cbuf,state->x,0,state->natoms);
1698 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1700 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1701 bUpdateDoLR,fr->f_twin,fcd,
1702 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1703 wallcycle_stop(wcycle,ewcUPDATE);
1705 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1706 /* are the small terms in the shake_vir here due
1707 * to numerical errors, or are they important
1708 * physically? I'm thinking they are just errors, but not completely sure.
1709 * For now, will call without actually constraining, constr=NULL*/
1710 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1711 state,fr->bMolPBC,graph,f,
1712 &top->idef,tmp_vir,force_vir,
1713 cr,nrnb,wcycle,upd,NULL,
1714 bInitStep,FALSE,bCalcVir,
1717 if (!bOK && !bFFscan)
1719 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1722 if (fr->bSepDVDL && fplog && do_log)
1724 fprintf(fplog,sepdvdlformat,"Constraint dV/dl",0.0,dvdl);
1726 enerd->term[F_DVDL_BONDED] += dvdl;
1730 /* Need to unshift here */
1731 unshift_self(graph,state->box,state->x);
1736 wallcycle_start(wcycle,ewcVSITECONSTR);
1739 shift_self(graph,state->box,state->x);
1741 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
1742 top->idef.iparams,top->idef.il,
1743 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
1747 unshift_self(graph,state->box,state->x);
1749 wallcycle_stop(wcycle,ewcVSITECONSTR);
1752 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1753 /* With Leap-Frog we can skip compute_globals at
1754 * non-communication steps, but we need to calculate
1755 * the kinetic energy one step before communication.
1757 if (bGStat || do_per_step(step+1,nstglobalcomm) ||
1760 if (ir->nstlist == -1 && bFirstIterate)
1762 gs.sig[eglsNABNSB] = nlh.nabnsb;
1764 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1765 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1767 bFirstIterate ? &gs : NULL,
1768 (step_rel % gs.nstms == 0) &&
1769 (multisim_nsteps<0 || (step_rel<multisim_nsteps)),
1771 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1773 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1774 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1775 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1776 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1777 | (bIterations && iterate.bIterate ? CGLO_ITERATE : 0)
1778 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1781 if (ir->nstlist == -1 && bFirstIterate)
1783 nlh.nabnsb = gs.set[eglsNABNSB];
1784 gs.set[eglsNABNSB] = 0;
1787 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1788 /* ############# END CALC EKIN AND PRESSURE ################# */
1790 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1791 the virial that should probably be addressed eventually. state->veta has better properies,
1792 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1793 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1796 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1797 trace(shake_vir),&tracevir))
1801 bFirstIterate = FALSE;
1804 /* only add constraint dvdl after constraints */
1805 enerd->term[F_DVDL_BONDED] += dvdl;
1806 if (!bVV || bRerunMD)
1808 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1809 sum_dhdl(enerd,state->lambda,ir->fepvals);
1811 update_box(fplog,step,ir,mdatoms,state,graph,f,
1812 ir->nstlist==-1 ? &nlh.scale_tot : NULL,pcoupl_mu,nrnb,wcycle,upd,bInitStep,FALSE);
1814 /* ################# END UPDATE STEP 2 ################# */
1815 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1817 /* The coordinates (x) were unshifted in update */
1818 if (bFFscan && (shellfc==NULL || bConverged))
1820 if (print_forcefield(fplog,enerd->term,mdatoms->homenr,
1822 &(top_global->mols),mdatoms->massT,pres))
1826 fprintf(stderr,"\n");
1832 /* We will not sum ekinh_old,
1833 * so signal that we still have to do it.
1835 bSumEkinhOld = TRUE;
1840 /* Only do GCT when the relaxation of shells (minimization) has converged,
1841 * otherwise we might be coupling to bogus energies.
1842 * In parallel we must always do this, because the other sims might
1846 /* Since this is called with the new coordinates state->x, I assume
1847 * we want the new box state->box too. / EL 20040121
1849 do_coupling(fplog,oenv,nfile,fnm,tcr,t,step,enerd->term,fr,
1851 mdatoms,&(top->idef),mu_aver,
1852 top_global->mols.nr,cr,
1853 state->box,total_vir,pres,
1854 mu_tot,state->x,f,bConverged);
1858 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1860 /* use the directly determined last velocity, not actually the averaged half steps */
1861 if (bTrotter && ir->eI==eiVV)
1863 enerd->term[F_EKIN] = last_ekin;
1865 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1869 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1873 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir,state,&MassQ);
1875 /* Check for excessively large energies */
1879 real etot_max = 1e200;
1881 real etot_max = 1e30;
1883 if (fabs(enerd->term[F_ETOT]) > etot_max)
1885 fprintf(stderr,"Energy too large (%g), giving up\n",
1886 enerd->term[F_ETOT]);
1889 /* ######### END PREPARING EDR OUTPUT ########### */
1891 /* Time for performance */
1892 if (((step % stepout) == 0) || bLastStep)
1894 runtime_upd_proc(runtime);
1900 gmx_bool do_dr,do_or;
1902 if (fplog && do_log && bDoExpanded)
1904 /* only needed if doing expanded ensemble */
1905 PrintFreeEnergyInfoToFile(fplog,ir->fepvals,ir->expandedvals,ir->bSimTemp?ir->simtempvals:NULL,
1906 &df_history,state->fep_state,ir->nstlog,step);
1908 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1912 upd_mdebin(mdebin,bDoDHDL, TRUE,
1913 t,mdatoms->tmass,enerd,state,
1914 ir->fepvals,ir->expandedvals,lastbox,
1915 shake_vir,force_vir,total_vir,pres,
1916 ekind,mu_tot,constr);
1920 upd_mdebin_step(mdebin);
1923 do_dr = do_per_step(step,ir->nstdisreout);
1924 do_or = do_per_step(step,ir->nstorireout);
1926 print_ebin(outf->fp_ene,do_ene,do_dr,do_or,do_log?fplog:NULL,
1928 eprNORMAL,bCompact,mdebin,fcd,groups,&(ir->opts));
1930 if (ir->ePull != epullNO)
1932 pull_print_output(ir->pull,step,t);
1935 if (do_per_step(step,ir->nstlog))
1937 if(fflush(fplog) != 0)
1939 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of disk space?");
1945 /* Have to do this part after outputting the logfile and the edr file */
1946 state->fep_state = lamnew;
1947 for (i=0;i<efptNR;i++)
1949 state_global->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
1952 /* Remaining runtime */
1953 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1957 fprintf(stderr,"\n");
1959 print_time(stderr,runtime,step,ir,cr);
1962 /* Replica exchange */
1964 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1965 do_per_step(step,repl_ex_nst))
1967 bExchanged = replica_exchange(fplog,cr,repl_ex,
1971 if (bExchanged && DOMAINDECOMP(cr))
1973 dd_partition_system(fplog,step,cr,TRUE,1,
1974 state_global,top_global,ir,
1975 state,&f,mdatoms,top,fr,
1976 vsite,shellfc,constr,
1983 bStartingFromCpt = FALSE;
1985 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1986 /* With all integrators, except VV, we need to retain the pressure
1987 * at the current step for coupling at the next step.
1989 if ((state->flags & (1<<estPRES_PREV)) &&
1991 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1993 /* Store the pressure in t_state for pressure coupling
1994 * at the next MD step.
1996 copy_mat(pres,state->pres_prev);
1999 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2001 if ( (membed!=NULL) && (!bLastStep) )
2003 rescale_membed(step_rel,membed,state_global->x);
2010 /* read next frame from input trajectory */
2011 bNotLastFrame = read_next_frame(oenv,status,&rerun_fr);
2016 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
2020 if (!bRerunMD || !rerun_fr.bStep)
2022 /* increase the MD step number */
2027 cycles = wallcycle_stop(wcycle,ewcSTEP);
2028 if (DOMAINDECOMP(cr) && wcycle)
2030 dd_cycles_add(cr->dd,cycles,ddCyclStep);
2033 if (bPMETuneRunning || bPMETuneTry)
2035 /* PME grid + cut-off optimization with GPUs or PME nodes */
2037 /* Count the total cycles over the last steps */
2038 cycles_pmes += cycles;
2040 /* We can only switch cut-off at NS steps */
2041 if (step % ir->nstlist == 0)
2043 /* PME grid + cut-off optimization with GPUs or PME nodes */
2046 if (DDMASTER(cr->dd))
2048 /* PME node load is too high, start tuning */
2049 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2051 dd_bcast(cr->dd,sizeof(gmx_bool),&bPMETuneRunning);
2053 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2055 bPMETuneTry = FALSE;
2058 if (bPMETuneRunning)
2060 /* init_step might not be a multiple of nstlist,
2061 * but the first cycle is always skipped anyhow.
2064 pme_load_balance(pme_loadbal,cr,
2065 (bVerbose && MASTER(cr)) ? stderr : NULL,
2067 ir,state,cycles_pmes,
2068 fr->ic,fr->nbv,&fr->pmedata,
2071 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2072 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2073 fr->rlist = fr->ic->rlist;
2074 fr->rlistlong = fr->ic->rlistlong;
2075 fr->rcoulomb = fr->ic->rcoulomb;
2076 fr->rvdw = fr->ic->rvdw;
2082 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2083 gs.set[eglsRESETCOUNTERS] != 0)
2085 /* Reset all the counters related to performance over the run */
2086 reset_all_counters(fplog,cr,step,&step_rel,ir,wcycle,nrnb,runtime,
2087 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2088 wcycle_set_reset_counters(wcycle,-1);
2089 /* Correct max_hours for the elapsed time */
2090 max_hours -= run_time/(60.0*60.0);
2091 bResetCountersHalfMaxH = FALSE;
2092 gs.set[eglsRESETCOUNTERS] = 0;
2096 /* End of main MD loop */
2100 runtime_end(runtime);
2102 if (bRerunMD && MASTER(cr))
2107 if (!(cr->duty & DUTY_PME))
2109 /* Tell the PME only node to finish */
2110 gmx_pme_send_finish(cr);
2115 if (ir->nstcalcenergy > 0 && !bRerunMD)
2117 print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
2118 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
2126 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2128 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)));
2129 fprintf(fplog,"Average number of atoms that crossed the half buffer length: %.1f\n\n",nlh.ab/nlh.nns);
2132 if (pme_loadbal != NULL)
2134 pme_loadbal_done(pme_loadbal,fplog);
2137 if (shellfc && fplog)
2139 fprintf(fplog,"Fraction of iterations that converged: %.2f %%\n",
2140 (nconverged*100.0)/step_rel);
2141 fprintf(fplog,"Average number of force evaluations per MD step: %.2f\n\n",
2145 if (repl_ex_nst > 0 && MASTER(cr))
2147 print_replica_exchange_statistics(fplog,repl_ex);
2150 runtime->nsteps_done = step_rel;