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"
72 #include "domdec_network.h"
78 #include "compute_io.h"
80 #include "checkpoint.h"
81 #include "mtop_util.h"
82 #include "sighandler.h"
85 #include "pme_loadbal.h"
88 #include "types/nlistheuristics.h"
89 #include "types/iteratedconstraints.h"
90 #include "nbnxn_cuda_data_mgmt.h"
100 #include "corewrap.h"
103 static void reset_all_counters(FILE *fplog,t_commrec *cr,
104 gmx_large_int_t step,
105 gmx_large_int_t *step_rel,t_inputrec *ir,
106 gmx_wallcycle_t wcycle,t_nrnb *nrnb,
107 gmx_runtime_t *runtime,
108 nbnxn_cuda_ptr_t cu_nbv)
110 char sbuf[STEPSTRSIZE];
112 /* Reset all the counters related to performance over the run */
113 md_print_warn(cr,fplog,"step %s: resetting all time and cycle counters\n",
114 gmx_step_str(step,sbuf));
118 nbnxn_cuda_reset_timings(cu_nbv);
121 wallcycle_stop(wcycle,ewcRUN);
122 wallcycle_reset_all(wcycle);
123 if (DOMAINDECOMP(cr))
125 reset_dd_statistics_counters(cr->dd);
128 ir->init_step += *step_rel;
129 ir->nsteps -= *step_rel;
131 wallcycle_start(wcycle,ewcRUN);
132 runtime_start(runtime);
133 print_date_and_time(fplog,cr->nodeid,"Restarted time",runtime);
136 double do_md(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
137 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
139 gmx_vsite_t *vsite,gmx_constr_t constr,
140 int stepout,t_inputrec *ir,
141 gmx_mtop_t *top_global,
143 t_state *state_global,
145 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
146 gmx_edsam_t ed,t_forcerec *fr,
147 int repl_ex_nst,int repl_ex_nex,int repl_ex_seed,gmx_membed_t membed,
148 real cpt_period,real max_hours,
149 const char *deviceOptions,
151 gmx_runtime_t *runtime)
154 gmx_large_int_t step,step_rel;
156 double t,t0,lam0[efptNR];
157 gmx_bool bGStatEveryStep,bGStat,bCalcVir,bCalcEner;
158 gmx_bool bNS,bNStList,bSimAnn,bStopCM,bRerunMD,bNotLastFrame=FALSE,
159 bFirstStep,bStateFromCP,bStateFromTPX,bInitStep,bLastStep,
160 bBornRadii,bStartingFromCpt;
161 gmx_bool bDoDHDL=FALSE,bDoFEP=FALSE,bDoExpanded=FALSE;
162 gmx_bool do_ene,do_log,do_verbose,bRerunWarnNoV=TRUE,
163 bForceUpdate=FALSE,bCPT;
165 gmx_bool bMasterState;
166 int force_flags,cglo_flags;
167 tensor force_vir,shake_vir,total_vir,tmp_vir,pres;
172 t_state *bufstate=NULL;
173 matrix *scale_tot,pcoupl_mu,M,ebox;
176 gmx_repl_ex_t repl_ex=NULL;
179 t_mdebin *mdebin=NULL;
180 df_history_t df_history;
185 gmx_enerdata_t *enerd;
187 gmx_global_stat_t gstat;
188 gmx_update_t upd=NULL;
191 gmx_rng_t mcrng=NULL;
193 gmx_groups_t *groups;
194 gmx_ekindata_t *ekind, *ekind_save;
195 gmx_shellfc_t shellfc;
196 int count,nconverged=0;
199 gmx_bool bIonize=FALSE;
200 gmx_bool bTCR=FALSE,bConverged=TRUE,bOK,bSumEkinhOld,bExchanged;
202 gmx_bool bResetCountersHalfMaxH=FALSE;
203 gmx_bool bVV,bIterations,bFirstIterate,bTemp,bPres,bTrotter;
206 atom_id *grpindex=NULL;
208 t_coupl_rec *tcr=NULL;
209 rvec *xcopy=NULL,*vcopy=NULL,*cbuf=NULL;
210 matrix boxcopy={{0}},lastbox;
212 real fom,oldfom,veta_save,pcurr,scalevir,tracevir;
220 real saved_conserved_quantity = 0;
225 char sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
226 int handled_stop_condition=gmx_stop_cond_none; /* compare to get_stop_condition*/
227 gmx_iterate_t iterate;
228 gmx_large_int_t multisim_nsteps=-1; /* number of steps to do before first multisim
229 simulation stops. If equal to zero, don't
230 communicate any more between multisims.*/
231 /* PME load balancing data for GPU kernels */
232 pme_load_balancing_t pme_loadbal=NULL;
234 gmx_bool bPMETuneTry=FALSE,bPMETuneRunning=FALSE;
239 "\n* WARNING * WARNING * WARNING * WARNING * WARNING * WARNING *\n"
240 "We have just committed the new CPU detection code in this branch,\n"
241 "and will commit new SSE/AVX kernels in a few days. However, this\n"
242 "means that currently only the NxN kernels are accelerated!\n"
243 "In the mean time, you might want to avoid production runs in 4.6.\n\n");
247 /* Temporary addition for FAHCORE checkpointing */
251 /* Check for special mdrun options */
252 bRerunMD = (Flags & MD_RERUN);
253 bIonize = (Flags & MD_IONIZE);
254 bFFscan = (Flags & MD_FFSCAN);
255 bAppend = (Flags & MD_APPENDFILES);
256 if (Flags & MD_RESETCOUNTERSHALFWAY)
260 /* Signal to reset the counters half the simulation steps. */
261 wcycle_set_reset_counters(wcycle,ir->nsteps/2);
263 /* Signal to reset the counters halfway the simulation time. */
264 bResetCountersHalfMaxH = (max_hours > 0);
267 /* md-vv uses averaged full step velocities for T-control
268 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
269 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
271 if (bVV) /* to store the initial velocities while computing virial */
273 snew(cbuf,top_global->natoms);
275 /* all the iteratative cases - only if there are constraints */
276 bIterations = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
277 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
281 /* Since we don't know if the frames read are related in any way,
282 * rebuild the neighborlist at every step.
285 ir->nstcalcenergy = 1;
289 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
291 nstglobalcomm = check_nstglobalcomm(fplog,cr,nstglobalcomm,ir);
292 bGStatEveryStep = (nstglobalcomm == 1);
294 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
297 "To reduce the energy communication with nstlist = -1\n"
298 "the neighbor list validity should not be checked at every step,\n"
299 "this means that exact integration is not guaranteed.\n"
300 "The neighbor list validity is checked after:\n"
301 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
302 "In most cases this will result in exact integration.\n"
303 "This reduces the energy communication by a factor of 2 to 3.\n"
304 "If you want less energy communication, set nstlist > 3.\n\n");
307 if (bRerunMD || bFFscan)
311 groups = &top_global->groups;
314 init_md(fplog,cr,ir,oenv,&t,&t0,state_global->lambda,
315 &(state_global->fep_state),lam0,
316 nrnb,top_global,&upd,
317 nfile,fnm,&outf,&mdebin,
318 force_vir,shake_vir,mu_tot,&bSimAnn,&vcm,state_global,Flags);
320 clear_mat(total_vir);
322 /* Energy terms and groups */
324 init_enerdata(top_global->groups.grps[egcENER].nr,ir->fepvals->n_lambda,
326 if (DOMAINDECOMP(cr))
332 snew(f,top_global->natoms);
335 /* lambda Monte carlo random number generator */
338 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
340 /* copy the state into df_history */
341 copy_df_history(&df_history,&state_global->dfhist);
343 /* Kinetic energy data */
345 init_ekindata(fplog,top_global,&(ir->opts),ekind);
346 /* needed for iteration of constraints */
348 init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
349 /* Copy the cos acceleration to the groups struct */
350 ekind->cosacc.cos_accel = ir->cos_accel;
352 gstat = global_stat_init(ir);
355 /* Check for polarizable models and flexible constraints */
356 shellfc = init_shell_flexcon(fplog,
357 top_global,n_flexible_constraints(constr),
358 (ir->bContinuation ||
359 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
360 NULL : state_global->x);
364 #ifdef GMX_THREAD_MPI
365 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
367 set_deform_reference_box(upd,
368 deform_init_init_step_tpx,
369 deform_init_box_tpx);
370 #ifdef GMX_THREAD_MPI
371 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
376 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
377 if ((io > 2000) && MASTER(cr))
379 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
383 if (DOMAINDECOMP(cr)) {
384 top = dd_init_local_top(top_global);
387 dd_init_local_state(cr->dd,state_global,state);
389 if (DDMASTER(cr->dd) && ir->nstfout) {
390 snew(f_global,state_global->natoms);
394 /* Initialize the particle decomposition and split the topology */
395 top = split_system(fplog,top_global,ir,cr);
397 pd_cg_range(cr,&fr->cg0,&fr->hcg);
398 pd_at_range(cr,&a0,&a1);
400 top = gmx_mtop_generate_local_top(top_global,ir);
403 a1 = top_global->natoms;
406 forcerec_set_excl_load(fr,top,cr);
408 state = partdec_init_local_state(cr,state_global);
411 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
414 set_vsite_top(vsite,top,mdatoms,cr);
417 if (ir->ePBC != epbcNONE && !fr->bMolPBC) {
418 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
422 make_local_shells(cr,mdatoms,shellfc);
425 init_bonded_thread_force_reduction(fr,&top->idef);
427 if (ir->pull && PAR(cr)) {
428 dd_make_local_pull_groups(NULL,ir->pull,mdatoms);
432 if (DOMAINDECOMP(cr))
434 /* Distribute the charge groups over the nodes from the master node */
435 dd_partition_system(fplog,ir->init_step,cr,TRUE,1,
436 state_global,top_global,ir,
437 state,&f,mdatoms,top,fr,
438 vsite,shellfc,constr,
443 update_mdatoms(mdatoms,state->lambda[efptMASS]);
445 if (opt2bSet("-cpi",nfile,fnm))
447 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr);
451 bStateFromCP = FALSE;
458 /* Update mdebin with energy history if appending to output files */
459 if ( Flags & MD_APPENDFILES )
461 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
465 /* We might have read an energy history from checkpoint,
466 * free the allocated memory and reset the counts.
468 done_energyhistory(&state_global->enerhist);
469 init_energyhistory(&state_global->enerhist);
472 /* Set the initial energy history in state by updating once */
473 update_energyhistory(&state_global->enerhist,mdebin);
476 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
478 /* Set the random state if we read a checkpoint file */
479 set_stochd_state(upd,state);
482 if (state->flags & (1<<estMC_RNG))
484 set_mc_state(mcrng,state);
487 /* Initialize constraints */
489 if (!DOMAINDECOMP(cr))
490 set_constraints(constr,top,ir,mdatoms,cr);
493 /* Check whether we have to GCT stuff */
494 bTCR = ftp2bSet(efGCT,nfile,fnm);
497 fprintf(stderr,"Will do General Coupling Theory!\n");
499 gnx = top_global->mols.nr;
501 for(i=0; (i<gnx); i++) {
508 /* We need to be sure replica exchange can only occur
509 * when the energies are current */
510 check_nst_param(fplog,cr,"nstcalcenergy",ir->nstcalcenergy,
511 "repl_ex_nst",&repl_ex_nst);
512 /* This check needs to happen before inter-simulation
513 * signals are initialized, too */
515 if (repl_ex_nst > 0 && MASTER(cr))
517 repl_ex = init_replica_exchange(fplog,cr->ms,state_global,ir,
518 repl_ex_nst,repl_ex_nex,repl_ex_seed);
521 /* PME tuning is only supported with GPUs or PME nodes and not with rerun */
522 if ((Flags & MD_TUNEPME) &&
523 EEL_PME(fr->eeltype) &&
524 fr->cutoff_scheme == ecutsVERLET &&
525 (fr->nbv->bUseGPU || !(cr->duty & DUTY_PME)) &&
528 pme_loadbal_init(&pme_loadbal,ir,state->box,fr->ic,fr->pmedata);
530 if (cr->duty & DUTY_PME)
532 /* Start tuning right away, as we can't measure the load */
533 bPMETuneRunning = TRUE;
537 /* Separate PME nodes, we can measure the PP/PME load balance */
542 if (!ir->bContinuation && !bRerunMD)
544 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
546 /* Set the velocities of frozen particles to zero */
547 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
551 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
561 /* Constrain the initial coordinates and velocities */
562 do_constrain_first(fplog,constr,ir,mdatoms,state,f,
563 graph,cr,nrnb,fr,top,shake_vir);
567 /* Construct the virtual sites for the initial configuration */
568 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
569 top->idef.iparams,top->idef.il,
570 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
576 /* set free energy calculation frequency as the minimum of nstdhdl, nstexpanded, and nstrepl_ex_nst*/
577 nstfep = ir->fepvals->nstdhdl;
578 if (ir->bExpanded && (nstfep > ir->expandedvals->nstexpanded))
580 nstfep = ir->expandedvals->nstexpanded;
582 if (repl_ex_nst > 0 && repl_ex_nst > nstfep)
584 nstfep = repl_ex_nst;
587 /* I'm assuming we need global communication the first time! MRS */
588 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
589 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM:0)
590 | (bVV ? CGLO_PRESSURE:0)
591 | (bVV ? CGLO_CONSTRAINT:0)
592 | (bRerunMD ? CGLO_RERUNMD:0)
593 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN:0));
595 bSumEkinhOld = FALSE;
596 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
597 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
598 constr,NULL,FALSE,state->box,
599 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,cglo_flags);
600 if (ir->eI == eiVVAK) {
601 /* a second call to get the half step temperature initialized as well */
602 /* we do the same call as above, but turn the pressure off -- internally to
603 compute_globals, this is recognized as a velocity verlet half-step
604 kinetic energy calculation. This minimized excess variables, but
605 perhaps loses some logic?*/
607 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
608 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
609 constr,NULL,FALSE,state->box,
610 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
611 cglo_flags &~ (CGLO_STOPCM | CGLO_PRESSURE));
614 /* Calculate the initial half step temperature, and save the ekinh_old */
615 if (!(Flags & MD_STARTFROMCPT))
617 for(i=0; (i<ir->opts.ngtc); i++)
619 copy_mat(ekind->tcstat[i].ekinh,ekind->tcstat[i].ekinh_old);
624 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
625 and there is no previous step */
628 /* if using an iterative algorithm, we need to create a working directory for the state. */
631 bufstate = init_bufstate(state);
635 snew(xcopy,state->natoms);
636 snew(vcopy,state->natoms);
637 copy_rvecn(state->x,xcopy,0,state->natoms);
638 copy_rvecn(state->v,vcopy,0,state->natoms);
639 copy_mat(state->box,boxcopy);
642 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
643 temperature control */
644 trotter_seq = init_npt_vars(ir,state,&MassQ,bTrotter);
648 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
651 "RMS relative constraint deviation after constraining: %.2e\n",
652 constr_rmsd(constr,FALSE));
654 if (EI_STATE_VELOCITY(ir->eI))
656 fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
660 fprintf(stderr,"starting md rerun '%s', reading coordinates from"
661 " input trajectory '%s'\n\n",
662 *(top_global->name),opt2fn("-rerun",nfile,fnm));
665 fprintf(stderr,"Calculated time to finish depends on nsteps from "
666 "run input file,\nwhich may not correspond to the time "
667 "needed to process input trajectory.\n\n");
673 fprintf(stderr,"starting mdrun '%s'\n",
674 *(top_global->name));
677 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
681 sprintf(tbuf,"%s","infinite");
683 if (ir->init_step > 0)
685 fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
686 gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
687 gmx_step_str(ir->init_step,sbuf2),
688 ir->init_step*ir->delta_t);
692 fprintf(stderr,"%s steps, %s ps.\n",
693 gmx_step_str(ir->nsteps,sbuf),tbuf);
699 /* Set and write start time */
700 runtime_start(runtime);
701 print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
702 wallcycle_start(wcycle,ewcRUN);
708 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
710 chkpt_ret=fcCheckPointParallel( cr->nodeid,
712 if ( chkpt_ret == 0 )
713 gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", 0 );
717 /***********************************************************
721 ************************************************************/
723 /* if rerunMD then read coordinates and velocities from input trajectory */
726 if (getenv("GMX_FORCE_UPDATE"))
734 bNotLastFrame = read_first_frame(oenv,&status,
735 opt2fn("-rerun",nfile,fnm),
736 &rerun_fr,TRX_NEED_X | TRX_READ_V);
737 if (rerun_fr.natoms != top_global->natoms)
740 "Number of atoms in trajectory (%d) does not match the "
741 "run input file (%d)\n",
742 rerun_fr.natoms,top_global->natoms);
744 if (ir->ePBC != epbcNONE)
748 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);
750 if (max_cutoff2(ir->ePBC,rerun_fr.box) < sqr(fr->rlistlong))
752 gmx_fatal(FARGS,"Rerun trajectory frame step %d time %f has too small box dimensions",rerun_fr.step,rerun_fr.time);
759 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
762 if (ir->ePBC != epbcNONE)
764 /* Set the shift vectors.
765 * Necessary here when have a static box different from the tpr box.
767 calc_shifts(rerun_fr.box,fr->shift_vec);
771 /* loop over MD steps or if rerunMD to end of input trajectory */
773 /* Skip the first Nose-Hoover integration when we get the state from tpx */
774 bStateFromTPX = !bStateFromCP;
775 bInitStep = bFirstStep && (bStateFromTPX || bVV);
776 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
778 bSumEkinhOld = FALSE;
781 init_global_signals(&gs,cr,ir,repl_ex_nst);
783 step = ir->init_step;
786 if (ir->nstlist == -1)
788 init_nlistheuristics(&nlh,bGStatEveryStep,step);
791 if (MULTISIM(cr) && (repl_ex_nst <=0 ))
793 /* check how many steps are left in other sims */
794 multisim_nsteps=get_multisim_nsteps(cr, ir->nsteps);
798 /* and stop now if we should */
799 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
800 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
801 while (!bLastStep || (bRerunMD && bNotLastFrame)) {
803 wallcycle_start(wcycle,ewcSTEP);
806 if (rerun_fr.bStep) {
807 step = rerun_fr.step;
808 step_rel = step - ir->init_step;
810 if (rerun_fr.bTime) {
820 bLastStep = (step_rel == ir->nsteps);
821 t = t0 + step*ir->delta_t;
824 if (ir->efep != efepNO || ir->bSimTemp)
826 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
827 requiring different logic. */
829 set_current_lambdas(step,ir->fepvals,bRerunMD,&rerun_fr,state_global,state,lam0);
830 bDoDHDL = do_per_step(step,ir->fepvals->nstdhdl);
831 bDoFEP = (do_per_step(step,nstfep) && (ir->efep != efepNO));
832 bDoExpanded = (do_per_step(step,ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
837 update_annealing_target_temp(&(ir->opts),t);
842 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
844 for(i=0; i<state_global->natoms; i++)
846 copy_rvec(rerun_fr.x[i],state_global->x[i]);
850 for(i=0; i<state_global->natoms; i++)
852 copy_rvec(rerun_fr.v[i],state_global->v[i]);
857 for(i=0; i<state_global->natoms; i++)
859 clear_rvec(state_global->v[i]);
863 fprintf(stderr,"\nWARNING: Some frames do not contain velocities.\n"
864 " Ekin, temperature and pressure are incorrect,\n"
865 " the virial will be incorrect when constraints are present.\n"
867 bRerunWarnNoV = FALSE;
871 copy_mat(rerun_fr.box,state_global->box);
872 copy_mat(state_global->box,state->box);
874 if (vsite && (Flags & MD_RERUN_VSITE))
876 if (DOMAINDECOMP(cr))
878 gmx_fatal(FARGS,"Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
882 /* Following is necessary because the graph may get out of sync
883 * with the coordinates if we only have every N'th coordinate set
885 mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
886 shift_self(graph,state->box,state->x);
888 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
889 top->idef.iparams,top->idef.il,
890 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
893 unshift_self(graph,state->box,state->x);
898 /* Stop Center of Mass motion */
899 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step,ir->nstcomm));
901 /* Copy back starting coordinates in case we're doing a forcefield scan */
904 for(ii=0; (ii<state->natoms); ii++)
906 copy_rvec(xcopy[ii],state->x[ii]);
907 copy_rvec(vcopy[ii],state->v[ii]);
909 copy_mat(boxcopy,state->box);
914 /* for rerun MD always do Neighbour Searching */
915 bNS = (bFirstStep || ir->nstlist != 0);
920 /* Determine whether or not to do Neighbour Searching and LR */
921 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
923 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
924 (ir->nstlist == -1 && nlh.nabnsb > 0));
926 if (bNS && ir->nstlist == -1)
928 set_nlistheuristics(&nlh,bFirstStep || bExchanged || bDoFEP, step);
932 /* check whether we should stop because another simulation has
936 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
937 (multisim_nsteps != ir->nsteps) )
944 "Stopping simulation %d because another one has finished\n",
948 gs.sig[eglsCHKPT] = 1;
953 /* < 0 means stop at next step, > 0 means stop at next NS step */
954 if ( (gs.set[eglsSTOPCOND] < 0 ) ||
955 ( (gs.set[eglsSTOPCOND] > 0 ) && ( bNS || ir->nstlist==0)) )
960 /* Determine whether or not to update the Born radii if doing GB */
961 bBornRadii=bFirstStep;
962 if (ir->implicit_solvent && (step % ir->nstgbradii==0))
967 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
968 do_verbose = bVerbose &&
969 (step % stepout == 0 || bFirstStep || bLastStep);
971 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
979 bMasterState = FALSE;
980 /* Correct the new box if it is too skewed */
981 if (DYNAMIC_BOX(*ir))
983 if (correct_box(fplog,step,state->box,graph))
988 if (DOMAINDECOMP(cr) && bMasterState)
990 dd_collect_state(cr->dd,state,state_global);
994 if (DOMAINDECOMP(cr))
996 /* Repartition the domain decomposition */
997 wallcycle_start(wcycle,ewcDOMDEC);
998 dd_partition_system(fplog,step,cr,
999 bMasterState,nstglobalcomm,
1000 state_global,top_global,ir,
1001 state,&f,mdatoms,top,fr,
1002 vsite,shellfc,constr,
1004 do_verbose && !bPMETuneRunning);
1005 wallcycle_stop(wcycle,ewcDOMDEC);
1006 /* If using an iterative integrator, reallocate space to match the decomposition */
1010 if (MASTER(cr) && do_log && !bFFscan)
1012 print_ebin_header(fplog,step,t,state->lambda[efptFEP]); /* can we improve the information printed here? */
1015 if (ir->efep != efepNO)
1017 update_mdatoms(mdatoms,state->lambda[efptMASS]);
1020 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1023 /* We need the kinetic energy at minus the half step for determining
1024 * the full step kinetic energy and possibly for T-coupling.*/
1025 /* This may not be quite working correctly yet . . . . */
1026 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1027 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1028 constr,NULL,FALSE,state->box,
1029 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1030 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1032 clear_mat(force_vir);
1034 /* Ionize the atoms if necessary */
1037 ionize(fplog,oenv,mdatoms,top_global,t,ir,state->x,state->v,
1038 mdatoms->start,mdatoms->start+mdatoms->homenr,state->box,cr);
1041 /* Update force field in ffscan program */
1044 if (update_forcefield(fplog,
1046 mdatoms->nr,state->x,state->box))
1054 /* We write a checkpoint at this MD step when:
1055 * either at an NS step when we signalled through gs,
1056 * or at the last step (but not when we do not want confout),
1057 * but never at the first step or with rerun.
1059 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1060 (bLastStep && (Flags & MD_CONFOUT))) &&
1061 step > ir->init_step && !bRerunMD);
1064 gs.set[eglsCHKPT] = 0;
1067 /* Determine the energy and pressure:
1068 * at nstcalcenergy steps and at energy output steps (set below).
1070 if (EI_VV(ir->eI) && (!bInitStep))
1072 /* for vv, the first half actually corresponds to the last step */
1073 bCalcEner = do_per_step(step-1,ir->nstcalcenergy);
1077 bCalcEner = do_per_step(step,ir->nstcalcenergy);
1079 bCalcVir = bCalcEner ||
1080 (ir->epc != epcNO && do_per_step(step,ir->nstpcouple));
1082 /* Do we need global communication ? */
1083 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1084 do_per_step(step,nstglobalcomm) ||
1085 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1087 do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
1089 if (do_ene || do_log)
1096 /* these CGLO_ options remain the same throughout the iteration */
1097 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1098 (bGStat ? CGLO_GSTAT : 0)
1101 force_flags = (GMX_FORCE_STATECHANGED |
1102 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1103 GMX_FORCE_ALLFORCES |
1104 (bNStList ? GMX_FORCE_DOLR : 0) |
1106 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1107 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1108 (bDoFEP ? GMX_FORCE_DHDL : 0)
1113 /* Now is the time to relax the shells */
1114 count=relax_shell_flexcon(fplog,cr,bVerbose,bFFscan ? step+1 : step,
1116 bStopCM,top,top_global,
1118 state,f,force_vir,mdatoms,
1119 nrnb,wcycle,graph,groups,
1120 shellfc,fr,bBornRadii,t,mu_tot,
1121 state->natoms,&bConverged,vsite,
1132 /* The coordinates (x) are shifted (to get whole molecules)
1134 * This is parallellized as well, and does communication too.
1135 * Check comments in sim_util.c
1137 do_force(fplog,cr,ir,step,nrnb,wcycle,top,top_global,groups,
1138 state->box,state->x,&state->hist,
1139 f,force_vir,mdatoms,enerd,fcd,
1140 state->lambda,graph,
1141 fr,vsite,mu_tot,t,outf->fp_field,ed,bBornRadii,
1142 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1147 mu_aver = calc_mu_aver(cr,state->x,mdatoms->chargeA,
1148 mu_tot,&top_global->mols,mdatoms,gnx,grpindex);
1151 if (bTCR && bFirstStep)
1153 tcr=init_coupling(fplog,nfile,fnm,cr,fr,mdatoms,&(top->idef));
1154 fprintf(fplog,"Done init_coupling\n");
1158 if (bVV && !bStartingFromCpt && !bRerunMD)
1159 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1161 if (ir->eI==eiVV && bInitStep)
1163 /* if using velocity verlet with full time step Ekin,
1164 * take the first half step only to compute the
1165 * virial for the first step. From there,
1166 * revert back to the initial coordinates
1167 * so that the input is actually the initial step.
1169 copy_rvecn(state->v,cbuf,0,state->natoms); /* should make this better for parallelizing? */
1171 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1172 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ1);
1175 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,
1176 f,fr->bTwinRange && bNStList,fr->f_twin,fcd,
1177 ekind,M,wcycle,upd,bInitStep,etrtVELOCITY1,
1178 cr,nrnb,constr,&top->idef);
1182 gmx_iterate_init(&iterate,bIterations && !bInitStep);
1184 /* for iterations, we save these vectors, as we will be self-consistently iterating
1187 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1189 /* save the state */
1190 if (bIterations && iterate.bIterate) {
1191 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1194 bFirstIterate = TRUE;
1195 while (bFirstIterate || (bIterations && iterate.bIterate))
1197 if (bIterations && iterate.bIterate)
1199 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1200 if (bFirstIterate && bTrotter)
1202 /* The first time through, we need a decent first estimate
1203 of veta(t+dt) to compute the constraints. Do
1204 this by computing the box volume part of the
1205 trotter integration at this time. Nothing else
1206 should be changed by this routine here. If
1207 !(first time), we start with the previous value
1210 veta_save = state->veta;
1211 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ0);
1212 vetanew = state->veta;
1213 state->veta = veta_save;
1218 if ( !bRerunMD || rerun_fr.bV || bForceUpdate) { /* Why is rerun_fr.bV here? Unclear. */
1221 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1222 state,fr->bMolPBC,graph,f,
1223 &top->idef,shake_vir,NULL,
1224 cr,nrnb,wcycle,upd,constr,
1225 bInitStep,TRUE,bCalcVir,vetanew);
1227 if (!bOK && !bFFscan)
1229 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1234 { /* Need to unshift here if a do_force has been
1235 called in the previous step */
1236 unshift_self(graph,state->box,state->x);
1240 /* if VV, compute the pressure and constraints */
1241 /* For VV2, we strictly only need this if using pressure
1242 * control, but we really would like to have accurate pressures
1244 * Think about ways around this in the future?
1245 * For now, keep this choice in comments.
1247 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1248 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1250 bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK));
1251 if (bCalcEner && ir->eI==eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1253 bSumEkinhOld = TRUE;
1255 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1256 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1257 constr,NULL,FALSE,state->box,
1258 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1261 | (bStopCM ? CGLO_STOPCM : 0)
1262 | (bTemp ? CGLO_TEMPERATURE:0)
1263 | (bPres ? CGLO_PRESSURE : 0)
1264 | (bPres ? CGLO_CONSTRAINT : 0)
1265 | ((bIterations && iterate.bIterate) ? CGLO_ITERATE : 0)
1266 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1269 /* explanation of above:
1270 a) We compute Ekin at the full time step
1271 if 1) we are using the AveVel Ekin, and it's not the
1272 initial step, or 2) if we are using AveEkin, but need the full
1273 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1274 b) If we are using EkinAveEkin for the kinetic energy for the temperture control, we still feed in
1275 EkinAveVel because it's needed for the pressure */
1277 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1282 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ2);
1289 /* We need the kinetic energy at minus the half step for determining
1290 * the full step kinetic energy and possibly for T-coupling.*/
1291 /* This may not be quite working correctly yet . . . . */
1292 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1293 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1294 constr,NULL,FALSE,state->box,
1295 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1296 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1300 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1305 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1306 state->veta,&vetanew))
1310 bFirstIterate = FALSE;
1313 if (bTrotter && !bInitStep) {
1314 enerd->term[F_DVDL_BONDED] += dvdl; /* only add after iterations */
1315 copy_mat(shake_vir,state->svir_prev);
1316 copy_mat(force_vir,state->fvir_prev);
1317 if (IR_NVT_TROTTER(ir) && ir->eI==eiVV) {
1318 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1319 enerd->term[F_TEMP] = sum_ekin(&(ir->opts),ekind,NULL,(ir->eI==eiVV),FALSE,FALSE);
1320 enerd->term[F_EKIN] = trace(ekind->ekin);
1323 /* if it's the initial step, we performed this first step just to get the constraint virial */
1324 if (bInitStep && ir->eI==eiVV) {
1325 copy_rvecn(cbuf,state->v,0,state->natoms);
1328 if (fr->bSepDVDL && fplog && do_log)
1330 fprintf(fplog,sepdvdlformat,"Constraint",0.0,dvdl);
1332 enerd->term[F_DVDL_BONDED] += dvdl;
1335 /* MRS -- now done iterating -- compute the conserved quantity */
1337 saved_conserved_quantity = compute_conserved_from_auxiliary(ir,state,&MassQ);
1340 last_ekin = enerd->term[F_EKIN];
1342 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1344 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1346 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1347 sum_dhdl(enerd,state->lambda,ir->fepvals);
1350 /* ######## END FIRST UPDATE STEP ############## */
1351 /* ######## If doing VV, we now have v(dt) ###### */
1353 /* perform extended ensemble sampling in lambda - we don't
1354 actually move to the new state before outputting
1355 statistics, but if performing simulated tempering, we
1356 do update the velocities and the tau_t. */
1358 lamnew = ExpandedEnsembleDynamics(fplog,ir,enerd,state,&MassQ,&df_history,step,mcrng,state->v,mdatoms);
1360 /* ################## START TRAJECTORY OUTPUT ################# */
1362 /* Now we have the energies and forces corresponding to the
1363 * coordinates at time t. We must output all of this before
1365 * for RerunMD t is read from input trajectory
1368 if (do_per_step(step,ir->nstxout)) { mdof_flags |= MDOF_X; }
1369 if (do_per_step(step,ir->nstvout)) { mdof_flags |= MDOF_V; }
1370 if (do_per_step(step,ir->nstfout)) { mdof_flags |= MDOF_F; }
1371 if (do_per_step(step,ir->nstxtcout)) { mdof_flags |= MDOF_XTC; }
1372 if (bCPT) { mdof_flags |= MDOF_CPT; };
1374 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1377 /* Enforce writing positions and velocities at end of run */
1378 mdof_flags |= (MDOF_X | MDOF_V);
1383 fcReportProgress( ir->nsteps, step );
1385 /* sync bCPT and fc record-keeping */
1386 if (bCPT && MASTER(cr))
1387 fcRequestCheckPoint();
1390 if (mdof_flags != 0)
1392 wallcycle_start(wcycle,ewcTRAJ);
1395 if (state->flags & (1<<estLD_RNG))
1397 get_stochd_state(upd,state);
1399 if (state->flags & (1<<estMC_RNG))
1401 get_mc_state(mcrng,state);
1407 state_global->ekinstate.bUpToDate = FALSE;
1411 update_ekinstate(&state_global->ekinstate,ekind);
1412 state_global->ekinstate.bUpToDate = TRUE;
1414 update_energyhistory(&state_global->enerhist,mdebin);
1415 if (ir->efep!=efepNO || ir->bSimTemp)
1417 state_global->fep_state = state->fep_state; /* MRS: seems kludgy. The code should be
1418 structured so this isn't necessary.
1419 Note this reassignment is only necessary
1420 for single threads.*/
1421 copy_df_history(&state_global->dfhist,&df_history);
1425 write_traj(fplog,cr,outf,mdof_flags,top_global,
1426 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
1433 if (bLastStep && step_rel == ir->nsteps &&
1434 (Flags & MD_CONFOUT) && MASTER(cr) &&
1435 !bRerunMD && !bFFscan)
1437 /* x and v have been collected in write_traj,
1438 * because a checkpoint file will always be written
1441 fprintf(stderr,"\nWriting final coordinates.\n");
1444 /* Make molecules whole only for confout writing */
1445 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
1447 write_sto_conf_mtop(ftp2fn(efSTO,nfile,fnm),
1448 *top_global->name,top_global,
1449 state_global->x,state_global->v,
1450 ir->ePBC,state->box);
1453 wallcycle_stop(wcycle,ewcTRAJ);
1456 /* kludge -- virial is lost with restart for NPT control. Must restart */
1457 if (bStartingFromCpt && bVV)
1459 copy_mat(state->svir_prev,shake_vir);
1460 copy_mat(state->fvir_prev,force_vir);
1462 /* ################## END TRAJECTORY OUTPUT ################ */
1464 /* Determine the wallclock run time up till now */
1465 run_time = gmx_gettime() - (double)runtime->real;
1467 /* Check whether everything is still allright */
1468 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1469 #ifdef GMX_THREAD_MPI
1474 /* this is just make gs.sig compatible with the hack
1475 of sending signals around by MPI_Reduce with together with
1477 if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
1478 gs.sig[eglsSTOPCOND]=1;
1479 if ( gmx_get_stop_condition() == gmx_stop_cond_next )
1480 gs.sig[eglsSTOPCOND]=-1;
1481 /* < 0 means stop at next step, > 0 means stop at next NS step */
1485 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1486 gmx_get_signal_name(),
1487 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1491 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1492 gmx_get_signal_name(),
1493 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1495 handled_stop_condition=(int)gmx_get_stop_condition();
1497 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1498 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1499 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1501 /* Signal to terminate the run */
1502 gs.sig[eglsSTOPCOND] = 1;
1505 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1507 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1510 if (bResetCountersHalfMaxH && MASTER(cr) &&
1511 run_time > max_hours*60.0*60.0*0.495)
1513 gs.sig[eglsRESETCOUNTERS] = 1;
1516 if (ir->nstlist == -1 && !bRerunMD)
1518 /* When bGStatEveryStep=FALSE, global_stat is only called
1519 * when we check the atom displacements, not at NS steps.
1520 * This means that also the bonded interaction count check is not
1521 * performed immediately after NS. Therefore a few MD steps could
1522 * be performed with missing interactions.
1523 * But wrong energies are never written to file,
1524 * since energies are only written after global_stat
1527 if (step >= nlh.step_nscheck)
1529 nlh.nabnsb = natoms_beyond_ns_buffer(ir,fr,&top->cgs,
1530 nlh.scale_tot,state->x);
1534 /* This is not necessarily true,
1535 * but step_nscheck is determined quite conservatively.
1541 /* In parallel we only have to check for checkpointing in steps
1542 * where we do global communication,
1543 * otherwise the other nodes don't know.
1545 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1548 run_time >= nchkpt*cpt_period*60.0)) &&
1549 gs.set[eglsCHKPT] == 0)
1551 gs.sig[eglsCHKPT] = 1;
1555 /* at the start of step, randomize the velocities */
1556 if (ETC_ANDERSEN(ir->etc) && EI_VV(ir->eI))
1558 gmx_bool bDoAndersenConstr;
1559 bDoAndersenConstr = (constr && update_randomize_velocities(ir,step,mdatoms,state,upd,&top->idef,constr));
1560 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1561 if (bDoAndersenConstr)
1563 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1564 state,fr->bMolPBC,graph,f,
1565 &top->idef,tmp_vir,NULL,
1566 cr,nrnb,wcycle,upd,constr,
1567 bInitStep,TRUE,bCalcVir,vetanew);
1573 gmx_iterate_init(&iterate,bIterations);
1576 /* for iterations, we save these vectors, as we will be redoing the calculations */
1577 if (bIterations && iterate.bIterate)
1579 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1581 bFirstIterate = TRUE;
1582 while (bFirstIterate || (bIterations && iterate.bIterate))
1584 /* We now restore these vectors to redo the calculation with improved extended variables */
1587 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1590 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1591 so scroll down for that logic */
1593 /* ######### START SECOND UPDATE STEP ################# */
1594 /* Box is changed in update() when we do pressure coupling,
1595 * but we should still use the old box for energy corrections and when
1596 * writing it to the energy file, so it matches the trajectory files for
1597 * the same timestep above. Make a copy in a separate array.
1599 copy_mat(state->box,lastbox);
1602 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1604 wallcycle_start(wcycle,ewcUPDATE);
1606 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1609 if (bIterations && iterate.bIterate)
1617 /* we use a new value of scalevir to converge the iterations faster */
1618 scalevir = tracevir/trace(shake_vir);
1620 msmul(shake_vir,scalevir,shake_vir);
1621 m_add(force_vir,shake_vir,total_vir);
1622 clear_mat(shake_vir);
1624 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ3);
1625 /* We can only do Berendsen coupling after we have summed
1626 * the kinetic energy or virial. Since the happens
1627 * in global_state after update, we should only do it at
1628 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1633 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1634 update_pcouple(fplog,step,ir,state,pcoupl_mu,M,wcycle,
1640 /* velocity half-step update */
1641 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1642 fr->bTwinRange && bNStList,fr->f_twin,fcd,
1643 ekind,M,wcycle,upd,FALSE,etrtVELOCITY2,
1644 cr,nrnb,constr,&top->idef);
1647 /* Above, initialize just copies ekinh into ekin,
1648 * it doesn't copy position (for VV),
1649 * and entire integrator for MD.
1654 copy_rvecn(state->x,cbuf,0,state->natoms);
1657 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1658 fr->bTwinRange && bNStList,fr->f_twin,fcd,
1659 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1660 wallcycle_stop(wcycle,ewcUPDATE);
1662 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,state,
1663 fr->bMolPBC,graph,f,
1664 &top->idef,shake_vir,force_vir,
1665 cr,nrnb,wcycle,upd,constr,
1666 bInitStep,FALSE,bCalcVir,state->veta);
1670 /* erase F_EKIN and F_TEMP here? */
1671 /* just compute the kinetic energy at the half step to perform a trotter step */
1672 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1673 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1674 constr,NULL,FALSE,lastbox,
1675 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1676 cglo_flags | CGLO_TEMPERATURE
1678 wallcycle_start(wcycle,ewcUPDATE);
1679 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ4);
1680 /* now we know the scaling, we can compute the positions again again */
1681 copy_rvecn(cbuf,state->x,0,state->natoms);
1683 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1684 fr->bTwinRange && bNStList,fr->f_twin,fcd,
1685 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1686 wallcycle_stop(wcycle,ewcUPDATE);
1688 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1689 /* are the small terms in the shake_vir here due
1690 * to numerical errors, or are they important
1691 * physically? I'm thinking they are just errors, but not completely sure.
1692 * For now, will call without actually constraining, constr=NULL*/
1693 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1694 state,fr->bMolPBC,graph,f,
1695 &top->idef,tmp_vir,force_vir,
1696 cr,nrnb,wcycle,upd,NULL,
1697 bInitStep,FALSE,bCalcVir,
1700 if (!bOK && !bFFscan)
1702 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1705 if (fr->bSepDVDL && fplog && do_log)
1707 fprintf(fplog,sepdvdlformat,"Constraint dV/dl",0.0,dvdl);
1709 enerd->term[F_DVDL_BONDED] += dvdl;
1713 /* Need to unshift here */
1714 unshift_self(graph,state->box,state->x);
1719 wallcycle_start(wcycle,ewcVSITECONSTR);
1722 shift_self(graph,state->box,state->x);
1724 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
1725 top->idef.iparams,top->idef.il,
1726 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
1730 unshift_self(graph,state->box,state->x);
1732 wallcycle_stop(wcycle,ewcVSITECONSTR);
1735 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1736 /* With Leap-Frog we can skip compute_globals at
1737 * non-communication steps, but we need to calculate
1738 * the kinetic energy one step before communication.
1740 if (bGStat || do_per_step(step+1,nstglobalcomm) ||
1743 if (ir->nstlist == -1 && bFirstIterate)
1745 gs.sig[eglsNABNSB] = nlh.nabnsb;
1747 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1748 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1750 bFirstIterate ? &gs : NULL,
1751 (step_rel % gs.nstms == 0) &&
1752 (multisim_nsteps<0 || (step_rel<multisim_nsteps)),
1754 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1756 | (!EI_VV(ir->eI) ? CGLO_ENERGY : 0)
1757 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1758 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1759 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1760 | (bIterations && iterate.bIterate ? CGLO_ITERATE : 0)
1761 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1764 if (ir->nstlist == -1 && bFirstIterate)
1766 nlh.nabnsb = gs.set[eglsNABNSB];
1767 gs.set[eglsNABNSB] = 0;
1770 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1771 /* ############# END CALC EKIN AND PRESSURE ################# */
1773 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1774 the virial that should probably be addressed eventually. state->veta has better properies,
1775 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1776 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1779 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1780 trace(shake_vir),&tracevir))
1784 bFirstIterate = FALSE;
1787 /* only add constraint dvdl after constraints */
1788 enerd->term[F_DVDL_BONDED] += dvdl;
1791 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1792 sum_dhdl(enerd,state->lambda,ir->fepvals);
1794 update_box(fplog,step,ir,mdatoms,state,graph,f,
1795 ir->nstlist==-1 ? &nlh.scale_tot : NULL,pcoupl_mu,nrnb,wcycle,upd,bInitStep,FALSE);
1797 /* ################# END UPDATE STEP 2 ################# */
1798 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1800 /* The coordinates (x) were unshifted in update */
1801 if (bFFscan && (shellfc==NULL || bConverged))
1803 if (print_forcefield(fplog,enerd->term,mdatoms->homenr,
1805 &(top_global->mols),mdatoms->massT,pres))
1809 fprintf(stderr,"\n");
1815 /* We will not sum ekinh_old,
1816 * so signal that we still have to do it.
1818 bSumEkinhOld = TRUE;
1823 /* Only do GCT when the relaxation of shells (minimization) has converged,
1824 * otherwise we might be coupling to bogus energies.
1825 * In parallel we must always do this, because the other sims might
1829 /* Since this is called with the new coordinates state->x, I assume
1830 * we want the new box state->box too. / EL 20040121
1832 do_coupling(fplog,oenv,nfile,fnm,tcr,t,step,enerd->term,fr,
1834 mdatoms,&(top->idef),mu_aver,
1835 top_global->mols.nr,cr,
1836 state->box,total_vir,pres,
1837 mu_tot,state->x,f,bConverged);
1841 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1843 /* use the directly determined last velocity, not actually the averaged half steps */
1844 if (bTrotter && ir->eI==eiVV)
1846 enerd->term[F_EKIN] = last_ekin;
1848 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1852 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1856 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir,state,&MassQ);
1858 /* Check for excessively large energies */
1862 real etot_max = 1e200;
1864 real etot_max = 1e30;
1866 if (fabs(enerd->term[F_ETOT]) > etot_max)
1868 fprintf(stderr,"Energy too large (%g), giving up\n",
1869 enerd->term[F_ETOT]);
1872 /* ######### END PREPARING EDR OUTPUT ########### */
1874 /* Time for performance */
1875 if (((step % stepout) == 0) || bLastStep)
1877 runtime_upd_proc(runtime);
1883 gmx_bool do_dr,do_or;
1885 if (fplog && do_log && bDoExpanded)
1887 /* only needed if doing expanded ensemble */
1888 PrintFreeEnergyInfoToFile(fplog,ir->fepvals,ir->expandedvals,ir->bSimTemp?ir->simtempvals:NULL,
1889 &df_history,state->fep_state,ir->nstlog,step);
1891 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1895 upd_mdebin(mdebin,bDoDHDL, TRUE,
1896 t,mdatoms->tmass,enerd,state,
1897 ir->fepvals,ir->expandedvals,lastbox,
1898 shake_vir,force_vir,total_vir,pres,
1899 ekind,mu_tot,constr);
1903 upd_mdebin_step(mdebin);
1906 do_dr = do_per_step(step,ir->nstdisreout);
1907 do_or = do_per_step(step,ir->nstorireout);
1909 print_ebin(outf->fp_ene,do_ene,do_dr,do_or,do_log?fplog:NULL,
1911 eprNORMAL,bCompact,mdebin,fcd,groups,&(ir->opts));
1913 if (ir->ePull != epullNO)
1915 pull_print_output(ir->pull,step,t);
1918 if (do_per_step(step,ir->nstlog))
1920 if(fflush(fplog) != 0)
1922 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of disk space?");
1928 /* Have to do this part after outputting the logfile and the edr file */
1929 state->fep_state = lamnew;
1930 for (i=0;i<efptNR;i++)
1932 state->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
1935 /* Remaining runtime */
1936 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1940 fprintf(stderr,"\n");
1942 print_time(stderr,runtime,step,ir,cr);
1945 /* Replica exchange */
1947 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1948 do_per_step(step,repl_ex_nst))
1950 bExchanged = replica_exchange(fplog,cr,repl_ex,
1954 if (bExchanged && DOMAINDECOMP(cr))
1956 dd_partition_system(fplog,step,cr,TRUE,1,
1957 state_global,top_global,ir,
1958 state,&f,mdatoms,top,fr,
1959 vsite,shellfc,constr,
1966 bStartingFromCpt = FALSE;
1968 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1969 /* With all integrators, except VV, we need to retain the pressure
1970 * at the current step for coupling at the next step.
1972 if ((state->flags & (1<<estPRES_PREV)) &&
1974 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1976 /* Store the pressure in t_state for pressure coupling
1977 * at the next MD step.
1979 copy_mat(pres,state->pres_prev);
1982 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1984 if ( (membed!=NULL) && (!bLastStep) )
1986 rescale_membed(step_rel,membed,state_global->x);
1993 /* read next frame from input trajectory */
1994 bNotLastFrame = read_next_frame(oenv,status,&rerun_fr);
1999 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
2003 if (!bRerunMD || !rerun_fr.bStep)
2005 /* increase the MD step number */
2010 cycles = wallcycle_stop(wcycle,ewcSTEP);
2011 if (DOMAINDECOMP(cr) && wcycle)
2013 dd_cycles_add(cr->dd,cycles,ddCyclStep);
2016 if (bPMETuneRunning || bPMETuneTry)
2018 /* PME grid + cut-off optimization with GPUs or PME nodes */
2020 /* Count the total cycles over the last steps */
2021 cycles_pmes += cycles;
2023 /* We can only switch cut-off at NS steps */
2024 if (step % ir->nstlist == 0)
2026 /* PME grid + cut-off optimization with GPUs or PME nodes */
2029 if (DDMASTER(cr->dd))
2031 /* PME node load is too high, start tuning */
2032 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2034 dd_bcast(cr->dd,sizeof(gmx_bool),&bPMETuneRunning);
2036 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2038 bPMETuneTry = FALSE;
2041 if (bPMETuneRunning)
2043 /* init_step might not be a multiple of nstlist,
2044 * but the first cycle is always skipped anyhow.
2047 pme_load_balance(pme_loadbal,cr,
2048 (bVerbose && MASTER(cr)) ? stderr : NULL,
2050 ir,state,cycles_pmes,
2051 fr->ic,fr->nbv,&fr->pmedata,
2054 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2061 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2062 gs.set[eglsRESETCOUNTERS] != 0)
2064 /* Reset all the counters related to performance over the run */
2065 reset_all_counters(fplog,cr,step,&step_rel,ir,wcycle,nrnb,runtime,
2066 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2067 wcycle_set_reset_counters(wcycle,-1);
2068 /* Correct max_hours for the elapsed time */
2069 max_hours -= run_time/(60.0*60.0);
2070 bResetCountersHalfMaxH = FALSE;
2071 gs.set[eglsRESETCOUNTERS] = 0;
2075 /* End of main MD loop */
2079 runtime_end(runtime);
2081 if (bRerunMD && MASTER(cr))
2086 if (!(cr->duty & DUTY_PME))
2088 /* Tell the PME only node to finish */
2089 gmx_pme_send_finish(cr);
2094 if (ir->nstcalcenergy > 0 && !bRerunMD)
2096 print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
2097 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
2105 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2107 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)));
2108 fprintf(fplog,"Average number of atoms that crossed the half buffer length: %.1f\n\n",nlh.ab/nlh.nns);
2111 if (pme_loadbal != NULL)
2113 pme_loadbal_done(pme_loadbal,fplog);
2116 if (shellfc && fplog)
2118 fprintf(fplog,"Fraction of iterations that converged: %.2f %%\n",
2119 (nconverged*100.0)/step_rel);
2120 fprintf(fplog,"Average number of force evaluations per MD step: %.2f\n\n",
2124 if (repl_ex_nst > 0 && MASTER(cr))
2126 print_replica_exchange_statistics(fplog,repl_ex);
2129 runtime->nsteps_done = step_rel;