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,bIterativeCase,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;
221 real saved_conserved_quantity = 0;
226 char sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
227 int handled_stop_condition=gmx_stop_cond_none; /* compare to get_stop_condition*/
228 gmx_iterate_t iterate;
229 gmx_large_int_t multisim_nsteps=-1; /* number of steps to do before first multisim
230 simulation stops. If equal to zero, don't
231 communicate any more between multisims.*/
232 /* PME load balancing data for GPU kernels */
233 pme_load_balancing_t pme_loadbal=NULL;
235 gmx_bool bPMETuneTry=FALSE,bPMETuneRunning=FALSE;
238 /* Temporary addition for FAHCORE checkpointing */
242 /* Check for special mdrun options */
243 bRerunMD = (Flags & MD_RERUN);
244 bIonize = (Flags & MD_IONIZE);
245 bFFscan = (Flags & MD_FFSCAN);
246 bAppend = (Flags & MD_APPENDFILES);
247 if (Flags & MD_RESETCOUNTERSHALFWAY)
251 /* Signal to reset the counters half the simulation steps. */
252 wcycle_set_reset_counters(wcycle,ir->nsteps/2);
254 /* Signal to reset the counters halfway the simulation time. */
255 bResetCountersHalfMaxH = (max_hours > 0);
258 /* md-vv uses averaged full step velocities for T-control
259 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
260 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
262 if (bVV) /* to store the initial velocities while computing virial */
264 snew(cbuf,top_global->natoms);
266 /* all the iteratative cases - only if there are constraints */
267 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
268 gmx_iterate_init(&iterate,FALSE); /* The default value of iterate->bIterationActive is set to
269 false in this step. The correct value, true or false,
270 is set at each step, as it depends on the frequency of temperature
271 and pressure control.*/
272 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
276 /* Since we don't know if the frames read are related in any way,
277 * rebuild the neighborlist at every step.
280 ir->nstcalcenergy = 1;
284 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
286 nstglobalcomm = check_nstglobalcomm(fplog,cr,nstglobalcomm,ir);
287 bGStatEveryStep = (nstglobalcomm == 1);
289 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
292 "To reduce the energy communication with nstlist = -1\n"
293 "the neighbor list validity should not be checked at every step,\n"
294 "this means that exact integration is not guaranteed.\n"
295 "The neighbor list validity is checked after:\n"
296 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
297 "In most cases this will result in exact integration.\n"
298 "This reduces the energy communication by a factor of 2 to 3.\n"
299 "If you want less energy communication, set nstlist > 3.\n\n");
302 if (bRerunMD || bFFscan)
306 groups = &top_global->groups;
309 init_md(fplog,cr,ir,oenv,&t,&t0,state_global->lambda,
310 &(state_global->fep_state),lam0,
311 nrnb,top_global,&upd,
312 nfile,fnm,&outf,&mdebin,
313 force_vir,shake_vir,mu_tot,&bSimAnn,&vcm,state_global,Flags);
315 clear_mat(total_vir);
317 /* Energy terms and groups */
319 init_enerdata(top_global->groups.grps[egcENER].nr,ir->fepvals->n_lambda,
321 if (DOMAINDECOMP(cr))
327 snew(f,top_global->natoms);
330 /* lambda Monte carlo random number generator */
333 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
335 /* copy the state into df_history */
336 copy_df_history(&df_history,&state_global->dfhist);
338 /* Kinetic energy data */
340 init_ekindata(fplog,top_global,&(ir->opts),ekind);
341 /* needed for iteration of constraints */
343 init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
344 /* Copy the cos acceleration to the groups struct */
345 ekind->cosacc.cos_accel = ir->cos_accel;
347 gstat = global_stat_init(ir);
350 /* Check for polarizable models and flexible constraints */
351 shellfc = init_shell_flexcon(fplog,
352 top_global,n_flexible_constraints(constr),
353 (ir->bContinuation ||
354 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
355 NULL : state_global->x);
359 #ifdef GMX_THREAD_MPI
360 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
362 set_deform_reference_box(upd,
363 deform_init_init_step_tpx,
364 deform_init_box_tpx);
365 #ifdef GMX_THREAD_MPI
366 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
371 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
372 if ((io > 2000) && MASTER(cr))
374 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
378 if (DOMAINDECOMP(cr)) {
379 top = dd_init_local_top(top_global);
382 dd_init_local_state(cr->dd,state_global,state);
384 if (DDMASTER(cr->dd) && ir->nstfout) {
385 snew(f_global,state_global->natoms);
389 /* Initialize the particle decomposition and split the topology */
390 top = split_system(fplog,top_global,ir,cr);
392 pd_cg_range(cr,&fr->cg0,&fr->hcg);
393 pd_at_range(cr,&a0,&a1);
395 top = gmx_mtop_generate_local_top(top_global,ir);
398 a1 = top_global->natoms;
401 forcerec_set_excl_load(fr,top,cr);
403 state = partdec_init_local_state(cr,state_global);
406 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
409 set_vsite_top(vsite,top,mdatoms,cr);
412 if (ir->ePBC != epbcNONE && !fr->bMolPBC) {
413 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
417 make_local_shells(cr,mdatoms,shellfc);
420 init_bonded_thread_force_reduction(fr,&top->idef);
422 if (ir->pull && PAR(cr)) {
423 dd_make_local_pull_groups(NULL,ir->pull,mdatoms);
427 if (DOMAINDECOMP(cr))
429 /* Distribute the charge groups over the nodes from the master node */
430 dd_partition_system(fplog,ir->init_step,cr,TRUE,1,
431 state_global,top_global,ir,
432 state,&f,mdatoms,top,fr,
433 vsite,shellfc,constr,
438 update_mdatoms(mdatoms,state->lambda[efptMASS]);
440 if (opt2bSet("-cpi",nfile,fnm))
442 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr);
446 bStateFromCP = FALSE;
453 /* Update mdebin with energy history if appending to output files */
454 if ( Flags & MD_APPENDFILES )
456 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
460 /* We might have read an energy history from checkpoint,
461 * free the allocated memory and reset the counts.
463 done_energyhistory(&state_global->enerhist);
464 init_energyhistory(&state_global->enerhist);
467 /* Set the initial energy history in state by updating once */
468 update_energyhistory(&state_global->enerhist,mdebin);
471 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
473 /* Set the random state if we read a checkpoint file */
474 set_stochd_state(upd,state);
477 if (state->flags & (1<<estMC_RNG))
479 set_mc_state(mcrng,state);
482 /* Initialize constraints */
484 if (!DOMAINDECOMP(cr))
485 set_constraints(constr,top,ir,mdatoms,cr);
488 /* Check whether we have to GCT stuff */
489 bTCR = ftp2bSet(efGCT,nfile,fnm);
492 fprintf(stderr,"Will do General Coupling Theory!\n");
494 gnx = top_global->mols.nr;
496 for(i=0; (i<gnx); i++) {
503 /* We need to be sure replica exchange can only occur
504 * when the energies are current */
505 check_nst_param(fplog,cr,"nstcalcenergy",ir->nstcalcenergy,
506 "repl_ex_nst",&repl_ex_nst);
507 /* This check needs to happen before inter-simulation
508 * signals are initialized, too */
510 if (repl_ex_nst > 0 && MASTER(cr))
512 repl_ex = init_replica_exchange(fplog,cr->ms,state_global,ir,
513 repl_ex_nst,repl_ex_nex,repl_ex_seed);
516 /* PME tuning is only supported with GPUs or PME nodes and not with rerun */
517 if ((Flags & MD_TUNEPME) &&
518 EEL_PME(fr->eeltype) &&
519 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
522 pme_loadbal_init(&pme_loadbal,ir,state->box,fr->ic,fr->pmedata);
524 if (cr->duty & DUTY_PME)
526 /* Start tuning right away, as we can't measure the load */
527 bPMETuneRunning = TRUE;
531 /* Separate PME nodes, we can measure the PP/PME load balance */
536 if (!ir->bContinuation && !bRerunMD)
538 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
540 /* Set the velocities of frozen particles to zero */
541 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
545 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
555 /* Constrain the initial coordinates and velocities */
556 do_constrain_first(fplog,constr,ir,mdatoms,state,f,
557 graph,cr,nrnb,fr,top,shake_vir);
561 /* Construct the virtual sites for the initial configuration */
562 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
563 top->idef.iparams,top->idef.il,
564 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
570 /* set free energy calculation frequency as the minimum of nstdhdl, nstexpanded, and nstrepl_ex_nst*/
571 nstfep = ir->fepvals->nstdhdl;
572 if (ir->bExpanded && (nstfep > ir->expandedvals->nstexpanded))
574 nstfep = ir->expandedvals->nstexpanded;
576 if (repl_ex_nst > 0 && nstfep > repl_ex_nst)
578 nstfep = repl_ex_nst;
581 /* I'm assuming we need global communication the first time! MRS */
582 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
583 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM:0)
584 | (bVV ? CGLO_PRESSURE:0)
585 | (bVV ? CGLO_CONSTRAINT:0)
586 | (bRerunMD ? CGLO_RERUNMD:0)
587 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN:0));
589 bSumEkinhOld = FALSE;
590 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
591 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
592 constr,NULL,FALSE,state->box,
593 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,cglo_flags);
594 if (ir->eI == eiVVAK) {
595 /* a second call to get the half step temperature initialized as well */
596 /* we do the same call as above, but turn the pressure off -- internally to
597 compute_globals, this is recognized as a velocity verlet half-step
598 kinetic energy calculation. This minimized excess variables, but
599 perhaps loses some logic?*/
601 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
602 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
603 constr,NULL,FALSE,state->box,
604 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
605 cglo_flags &~ (CGLO_STOPCM | CGLO_PRESSURE));
608 /* Calculate the initial half step temperature, and save the ekinh_old */
609 if (!(Flags & MD_STARTFROMCPT))
611 for(i=0; (i<ir->opts.ngtc); i++)
613 copy_mat(ekind->tcstat[i].ekinh,ekind->tcstat[i].ekinh_old);
618 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
619 and there is no previous step */
622 /* if using an iterative algorithm, we need to create a working directory for the state. */
625 bufstate = init_bufstate(state);
629 snew(xcopy,state->natoms);
630 snew(vcopy,state->natoms);
631 copy_rvecn(state->x,xcopy,0,state->natoms);
632 copy_rvecn(state->v,vcopy,0,state->natoms);
633 copy_mat(state->box,boxcopy);
636 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
637 temperature control */
638 trotter_seq = init_npt_vars(ir,state,&MassQ,bTrotter);
642 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
645 "RMS relative constraint deviation after constraining: %.2e\n",
646 constr_rmsd(constr,FALSE));
648 if (EI_STATE_VELOCITY(ir->eI))
650 fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
654 fprintf(stderr,"starting md rerun '%s', reading coordinates from"
655 " input trajectory '%s'\n\n",
656 *(top_global->name),opt2fn("-rerun",nfile,fnm));
659 fprintf(stderr,"Calculated time to finish depends on nsteps from "
660 "run input file,\nwhich may not correspond to the time "
661 "needed to process input trajectory.\n\n");
667 fprintf(stderr,"starting mdrun '%s'\n",
668 *(top_global->name));
671 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
675 sprintf(tbuf,"%s","infinite");
677 if (ir->init_step > 0)
679 fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
680 gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
681 gmx_step_str(ir->init_step,sbuf2),
682 ir->init_step*ir->delta_t);
686 fprintf(stderr,"%s steps, %s ps.\n",
687 gmx_step_str(ir->nsteps,sbuf),tbuf);
693 /* Set and write start time */
694 runtime_start(runtime);
695 print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
696 wallcycle_start(wcycle,ewcRUN);
702 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
704 chkpt_ret=fcCheckPointParallel( cr->nodeid,
706 if ( chkpt_ret == 0 )
707 gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", 0 );
711 /***********************************************************
715 ************************************************************/
717 /* if rerunMD then read coordinates and velocities from input trajectory */
720 if (getenv("GMX_FORCE_UPDATE"))
728 bNotLastFrame = read_first_frame(oenv,&status,
729 opt2fn("-rerun",nfile,fnm),
730 &rerun_fr,TRX_NEED_X | TRX_READ_V);
731 if (rerun_fr.natoms != top_global->natoms)
734 "Number of atoms in trajectory (%d) does not match the "
735 "run input file (%d)\n",
736 rerun_fr.natoms,top_global->natoms);
738 if (ir->ePBC != epbcNONE)
742 gmx_fatal(FARGS,"Rerun trajectory frame step %d time %f does not contain a box, while pbc is used",rerun_fr.step,rerun_fr.time);
744 if (max_cutoff2(ir->ePBC,rerun_fr.box) < sqr(fr->rlistlong))
746 gmx_fatal(FARGS,"Rerun trajectory frame step %d time %f has too small box dimensions",rerun_fr.step,rerun_fr.time);
753 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
756 if (ir->ePBC != epbcNONE)
758 /* Set the shift vectors.
759 * Necessary here when have a static box different from the tpr box.
761 calc_shifts(rerun_fr.box,fr->shift_vec);
765 /* loop over MD steps or if rerunMD to end of input trajectory */
767 /* Skip the first Nose-Hoover integration when we get the state from tpx */
768 bStateFromTPX = !bStateFromCP;
769 bInitStep = bFirstStep && (bStateFromTPX || bVV);
770 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
772 bSumEkinhOld = FALSE;
775 init_global_signals(&gs,cr,ir,repl_ex_nst);
777 step = ir->init_step;
780 if (ir->nstlist == -1)
782 init_nlistheuristics(&nlh,bGStatEveryStep,step);
785 if (MULTISIM(cr) && (repl_ex_nst <=0 ))
787 /* check how many steps are left in other sims */
788 multisim_nsteps=get_multisim_nsteps(cr, ir->nsteps);
792 /* and stop now if we should */
793 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
794 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
795 while (!bLastStep || (bRerunMD && bNotLastFrame)) {
797 wallcycle_start(wcycle,ewcSTEP);
800 if (rerun_fr.bStep) {
801 step = rerun_fr.step;
802 step_rel = step - ir->init_step;
804 if (rerun_fr.bTime) {
814 bLastStep = (step_rel == ir->nsteps);
815 t = t0 + step*ir->delta_t;
818 if (ir->efep != efepNO || ir->bSimTemp)
820 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
821 requiring different logic. */
823 set_current_lambdas(step,ir->fepvals,bRerunMD,&rerun_fr,state_global,state,lam0);
824 bDoDHDL = do_per_step(step,ir->fepvals->nstdhdl);
825 bDoFEP = (do_per_step(step,nstfep) && (ir->efep != efepNO));
826 bDoExpanded = (do_per_step(step,ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
831 update_annealing_target_temp(&(ir->opts),t);
836 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
838 for(i=0; i<state_global->natoms; i++)
840 copy_rvec(rerun_fr.x[i],state_global->x[i]);
844 for(i=0; i<state_global->natoms; i++)
846 copy_rvec(rerun_fr.v[i],state_global->v[i]);
851 for(i=0; i<state_global->natoms; i++)
853 clear_rvec(state_global->v[i]);
857 fprintf(stderr,"\nWARNING: Some frames do not contain velocities.\n"
858 " Ekin, temperature and pressure are incorrect,\n"
859 " the virial will be incorrect when constraints are present.\n"
861 bRerunWarnNoV = FALSE;
865 copy_mat(rerun_fr.box,state_global->box);
866 copy_mat(state_global->box,state->box);
868 if (vsite && (Flags & MD_RERUN_VSITE))
870 if (DOMAINDECOMP(cr))
872 gmx_fatal(FARGS,"Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
876 /* Following is necessary because the graph may get out of sync
877 * with the coordinates if we only have every N'th coordinate set
879 mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
880 shift_self(graph,state->box,state->x);
882 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
883 top->idef.iparams,top->idef.il,
884 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
887 unshift_self(graph,state->box,state->x);
892 /* Stop Center of Mass motion */
893 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step,ir->nstcomm));
895 /* Copy back starting coordinates in case we're doing a forcefield scan */
898 for(ii=0; (ii<state->natoms); ii++)
900 copy_rvec(xcopy[ii],state->x[ii]);
901 copy_rvec(vcopy[ii],state->v[ii]);
903 copy_mat(boxcopy,state->box);
908 /* for rerun MD always do Neighbour Searching */
909 bNS = (bFirstStep || ir->nstlist != 0);
914 /* Determine whether or not to do Neighbour Searching and LR */
915 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
917 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
918 (ir->nstlist == -1 && nlh.nabnsb > 0));
920 if (bNS && ir->nstlist == -1)
922 set_nlistheuristics(&nlh,bFirstStep || bExchanged || bDoFEP, step);
926 /* check whether we should stop because another simulation has
930 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
931 (multisim_nsteps != ir->nsteps) )
938 "Stopping simulation %d because another one has finished\n",
942 gs.sig[eglsCHKPT] = 1;
947 /* < 0 means stop at next step, > 0 means stop at next NS step */
948 if ( (gs.set[eglsSTOPCOND] < 0 ) ||
949 ( (gs.set[eglsSTOPCOND] > 0 ) && ( bNS || ir->nstlist==0)) )
954 /* Determine whether or not to update the Born radii if doing GB */
955 bBornRadii=bFirstStep;
956 if (ir->implicit_solvent && (step % ir->nstgbradii==0))
961 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
962 do_verbose = bVerbose &&
963 (step % stepout == 0 || bFirstStep || bLastStep);
965 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
973 bMasterState = FALSE;
974 /* Correct the new box if it is too skewed */
975 if (DYNAMIC_BOX(*ir))
977 if (correct_box(fplog,step,state->box,graph))
982 if (DOMAINDECOMP(cr) && bMasterState)
984 dd_collect_state(cr->dd,state,state_global);
988 if (DOMAINDECOMP(cr))
990 /* Repartition the domain decomposition */
991 wallcycle_start(wcycle,ewcDOMDEC);
992 dd_partition_system(fplog,step,cr,
993 bMasterState,nstglobalcomm,
994 state_global,top_global,ir,
995 state,&f,mdatoms,top,fr,
996 vsite,shellfc,constr,
998 do_verbose && !bPMETuneRunning);
999 wallcycle_stop(wcycle,ewcDOMDEC);
1000 /* If using an iterative integrator, reallocate space to match the decomposition */
1004 if (MASTER(cr) && do_log && !bFFscan)
1006 print_ebin_header(fplog,step,t,state->lambda[efptFEP]); /* can we improve the information printed here? */
1009 if (ir->efep != efepNO)
1011 update_mdatoms(mdatoms,state->lambda[efptMASS]);
1014 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1017 /* We need the kinetic energy at minus the half step for determining
1018 * the full step kinetic energy and possibly for T-coupling.*/
1019 /* This may not be quite working correctly yet . . . . */
1020 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1021 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1022 constr,NULL,FALSE,state->box,
1023 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1024 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1026 clear_mat(force_vir);
1028 /* Ionize the atoms if necessary */
1031 ionize(fplog,oenv,mdatoms,top_global,t,ir,state->x,state->v,
1032 mdatoms->start,mdatoms->start+mdatoms->homenr,state->box,cr);
1035 /* Update force field in ffscan program */
1038 if (update_forcefield(fplog,
1040 mdatoms->nr,state->x,state->box))
1048 /* We write a checkpoint at this MD step when:
1049 * either at an NS step when we signalled through gs,
1050 * or at the last step (but not when we do not want confout),
1051 * but never at the first step or with rerun.
1053 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1054 (bLastStep && (Flags & MD_CONFOUT))) &&
1055 step > ir->init_step && !bRerunMD);
1058 gs.set[eglsCHKPT] = 0;
1061 /* Determine the energy and pressure:
1062 * at nstcalcenergy steps and at energy output steps (set below).
1064 if (EI_VV(ir->eI) && (!bInitStep))
1066 /* for vv, the first half of the integration actually corresponds
1067 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1068 but the virial needs to be calculated on both the current step and the 'next' step. Future
1069 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1071 bCalcEner = do_per_step(step-1,ir->nstcalcenergy);
1072 bCalcVir = bCalcEner ||
1073 (ir->epc != epcNO && (do_per_step(step,ir->nstpcouple) || do_per_step(step-1,ir->nstpcouple)));
1077 bCalcEner = do_per_step(step,ir->nstcalcenergy);
1078 bCalcVir = bCalcEner ||
1079 (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 |
1105 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1106 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1107 (bDoFEP ? GMX_FORCE_DHDL : 0)
1112 if(do_per_step(step,ir->nstcalclr))
1114 force_flags |= GMX_FORCE_DO_LR;
1120 /* Now is the time to relax the shells */
1121 count=relax_shell_flexcon(fplog,cr,bVerbose,bFFscan ? step+1 : step,
1123 bStopCM,top,top_global,
1125 state,f,force_vir,mdatoms,
1126 nrnb,wcycle,graph,groups,
1127 shellfc,fr,bBornRadii,t,mu_tot,
1128 state->natoms,&bConverged,vsite,
1139 /* The coordinates (x) are shifted (to get whole molecules)
1141 * This is parallellized as well, and does communication too.
1142 * Check comments in sim_util.c
1144 do_force(fplog,cr,ir,step,nrnb,wcycle,top,top_global,groups,
1145 state->box,state->x,&state->hist,
1146 f,force_vir,mdatoms,enerd,fcd,
1147 state->lambda,graph,
1148 fr,vsite,mu_tot,t,outf->fp_field,ed,bBornRadii,
1149 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1154 mu_aver = calc_mu_aver(cr,state->x,mdatoms->chargeA,
1155 mu_tot,&top_global->mols,mdatoms,gnx,grpindex);
1158 if (bTCR && bFirstStep)
1160 tcr=init_coupling(fplog,nfile,fnm,cr,fr,mdatoms,&(top->idef));
1161 fprintf(fplog,"Done init_coupling\n");
1165 if (bVV && !bStartingFromCpt && !bRerunMD)
1166 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1168 if (ir->eI==eiVV && bInitStep)
1170 /* if using velocity verlet with full time step Ekin,
1171 * take the first half step only to compute the
1172 * virial for the first step. From there,
1173 * revert back to the initial coordinates
1174 * so that the input is actually the initial step.
1176 copy_rvecn(state->v,cbuf,0,state->natoms); /* should make this better for parallelizing? */
1178 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1179 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ1);
1182 /* If we are using twin-range interactions where the long-range component
1183 * is only evaluated every nstcalclr>1 steps, we should do a special update
1184 * step to combine the long-range forces on these steps.
1185 * For nstcalclr=1 this is not done, since the forces would have been added
1186 * directly to the short-range forces already.
1188 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1190 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,
1191 f,bUpdateDoLR,fr->f_twin,fcd,
1192 ekind,M,wcycle,upd,bInitStep,etrtVELOCITY1,
1193 cr,nrnb,constr,&top->idef);
1195 if (bIterativeCase && do_per_step(step-1,ir->nstpcouple) && !bInitStep)
1197 gmx_iterate_init(&iterate,TRUE);
1199 /* for iterations, we save these vectors, as we will be self-consistently iterating
1202 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1204 /* save the state */
1205 if (iterate.bIterationActive) {
1206 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1209 bFirstIterate = TRUE;
1210 while (bFirstIterate || iterate.bIterationActive)
1212 if (iterate.bIterationActive)
1214 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1215 if (bFirstIterate && bTrotter)
1217 /* The first time through, we need a decent first estimate
1218 of veta(t+dt) to compute the constraints. Do
1219 this by computing the box volume part of the
1220 trotter integration at this time. Nothing else
1221 should be changed by this routine here. If
1222 !(first time), we start with the previous value
1225 veta_save = state->veta;
1226 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ0);
1227 vetanew = state->veta;
1228 state->veta = veta_save;
1233 if ( !bRerunMD || rerun_fr.bV || bForceUpdate) { /* Why is rerun_fr.bV here? Unclear. */
1236 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1237 state,fr->bMolPBC,graph,f,
1238 &top->idef,shake_vir,NULL,
1239 cr,nrnb,wcycle,upd,constr,
1240 bInitStep,TRUE,bCalcVir,vetanew);
1242 if (!bOK && !bFFscan)
1244 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1249 { /* Need to unshift here if a do_force has been
1250 called in the previous step */
1251 unshift_self(graph,state->box,state->x);
1254 /* if VV, compute the pressure and constraints */
1255 /* For VV2, we strictly only need this if using pressure
1256 * control, but we really would like to have accurate pressures
1258 * Think about ways around this in the future?
1259 * For now, keep this choice in comments.
1261 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1262 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1264 bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK));
1265 if (bCalcEner && ir->eI==eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1267 bSumEkinhOld = TRUE;
1269 /* for vv, the first half of the integration actually corresponds to the previous step.
1270 So we need information from the last step in the first half of the integration */
1271 if (bGStat || do_per_step(step-1,nstglobalcomm)) {
1272 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1273 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1274 constr,NULL,FALSE,state->box,
1275 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1278 | (bTemp ? CGLO_TEMPERATURE:0)
1279 | (bPres ? CGLO_PRESSURE : 0)
1280 | (bPres ? CGLO_CONSTRAINT : 0)
1281 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1282 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1285 /* explanation of above:
1286 a) We compute Ekin at the full time step
1287 if 1) we are using the AveVel Ekin, and it's not the
1288 initial step, or 2) if we are using AveEkin, but need the full
1289 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1290 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1291 EkinAveVel because it's needed for the pressure */
1293 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1298 m_add(force_vir,shake_vir,total_vir); /* we need the un-dispersion corrected total vir here */
1299 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ2);
1306 /* We need the kinetic energy at minus the half step for determining
1307 * the full step kinetic energy and possibly for T-coupling.*/
1308 /* This may not be quite working correctly yet . . . . */
1309 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1310 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1311 constr,NULL,FALSE,state->box,
1312 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1313 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1318 if (iterate.bIterationActive &&
1319 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1320 state->veta,&vetanew))
1324 bFirstIterate = FALSE;
1327 if (bTrotter && !bInitStep) {
1328 enerd->term[F_DVDL_BONDED] += dvdl; /* only add after iterations */
1329 copy_mat(shake_vir,state->svir_prev);
1330 copy_mat(force_vir,state->fvir_prev);
1331 if (IR_NVT_TROTTER(ir) && ir->eI==eiVV) {
1332 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1333 enerd->term[F_TEMP] = sum_ekin(&(ir->opts),ekind,NULL,(ir->eI==eiVV),FALSE,FALSE);
1334 enerd->term[F_EKIN] = trace(ekind->ekin);
1337 /* if it's the initial step, we performed this first step just to get the constraint virial */
1338 if (bInitStep && ir->eI==eiVV) {
1339 copy_rvecn(cbuf,state->v,0,state->natoms);
1342 if (fr->bSepDVDL && fplog && do_log)
1344 fprintf(fplog,sepdvdlformat,"Constraint",0.0,dvdl);
1346 enerd->term[F_DVDL_BONDED] += dvdl;
1349 /* MRS -- now done iterating -- compute the conserved quantity */
1351 saved_conserved_quantity = compute_conserved_from_auxiliary(ir,state,&MassQ);
1354 last_ekin = enerd->term[F_EKIN];
1356 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1358 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1360 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1363 sum_dhdl(enerd,state->lambda,ir->fepvals);
1367 /* ######## END FIRST UPDATE STEP ############## */
1368 /* ######## If doing VV, we now have v(dt) ###### */
1370 /* perform extended ensemble sampling in lambda - we don't
1371 actually move to the new state before outputting
1372 statistics, but if performing simulated tempering, we
1373 do update the velocities and the tau_t. */
1375 lamnew = ExpandedEnsembleDynamics(fplog,ir,enerd,state,&MassQ,&df_history,step,mcrng,state->v,mdatoms);
1377 /* ################## START TRAJECTORY OUTPUT ################# */
1379 /* Now we have the energies and forces corresponding to the
1380 * coordinates at time t. We must output all of this before
1382 * for RerunMD t is read from input trajectory
1385 if (do_per_step(step,ir->nstxout)) { mdof_flags |= MDOF_X; }
1386 if (do_per_step(step,ir->nstvout)) { mdof_flags |= MDOF_V; }
1387 if (do_per_step(step,ir->nstfout)) { mdof_flags |= MDOF_F; }
1388 if (do_per_step(step,ir->nstxtcout)) { mdof_flags |= MDOF_XTC; }
1389 if (bCPT) { mdof_flags |= MDOF_CPT; };
1391 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1394 /* Enforce writing positions and velocities at end of run */
1395 mdof_flags |= (MDOF_X | MDOF_V);
1400 fcReportProgress( ir->nsteps, step );
1402 /* sync bCPT and fc record-keeping */
1403 if (bCPT && MASTER(cr))
1404 fcRequestCheckPoint();
1407 if (mdof_flags != 0)
1409 wallcycle_start(wcycle,ewcTRAJ);
1412 if (state->flags & (1<<estLD_RNG))
1414 get_stochd_state(upd,state);
1416 if (state->flags & (1<<estMC_RNG))
1418 get_mc_state(mcrng,state);
1424 state_global->ekinstate.bUpToDate = FALSE;
1428 update_ekinstate(&state_global->ekinstate,ekind);
1429 state_global->ekinstate.bUpToDate = TRUE;
1431 update_energyhistory(&state_global->enerhist,mdebin);
1432 if (ir->efep!=efepNO || ir->bSimTemp)
1434 state_global->fep_state = state->fep_state; /* MRS: seems kludgy. The code should be
1435 structured so this isn't necessary.
1436 Note this reassignment is only necessary
1437 for single threads.*/
1438 copy_df_history(&state_global->dfhist,&df_history);
1442 write_traj(fplog,cr,outf,mdof_flags,top_global,
1443 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
1450 if (bLastStep && step_rel == ir->nsteps &&
1451 (Flags & MD_CONFOUT) && MASTER(cr) &&
1452 !bRerunMD && !bFFscan)
1454 /* x and v have been collected in write_traj,
1455 * because a checkpoint file will always be written
1458 fprintf(stderr,"\nWriting final coordinates.\n");
1461 /* Make molecules whole only for confout writing */
1462 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
1464 write_sto_conf_mtop(ftp2fn(efSTO,nfile,fnm),
1465 *top_global->name,top_global,
1466 state_global->x,state_global->v,
1467 ir->ePBC,state->box);
1470 wallcycle_stop(wcycle,ewcTRAJ);
1473 /* kludge -- virial is lost with restart for NPT control. Must restart */
1474 if (bStartingFromCpt && bVV)
1476 copy_mat(state->svir_prev,shake_vir);
1477 copy_mat(state->fvir_prev,force_vir);
1479 /* ################## END TRAJECTORY OUTPUT ################ */
1481 /* Determine the wallclock run time up till now */
1482 run_time = gmx_gettime() - (double)runtime->real;
1484 /* Check whether everything is still allright */
1485 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1486 #ifdef GMX_THREAD_MPI
1491 /* this is just make gs.sig compatible with the hack
1492 of sending signals around by MPI_Reduce with together with
1494 if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
1495 gs.sig[eglsSTOPCOND]=1;
1496 if ( gmx_get_stop_condition() == gmx_stop_cond_next )
1497 gs.sig[eglsSTOPCOND]=-1;
1498 /* < 0 means stop at next step, > 0 means stop at next NS step */
1502 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1503 gmx_get_signal_name(),
1504 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1508 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1509 gmx_get_signal_name(),
1510 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1512 handled_stop_condition=(int)gmx_get_stop_condition();
1514 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1515 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1516 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1518 /* Signal to terminate the run */
1519 gs.sig[eglsSTOPCOND] = 1;
1522 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1524 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1527 if (bResetCountersHalfMaxH && MASTER(cr) &&
1528 run_time > max_hours*60.0*60.0*0.495)
1530 gs.sig[eglsRESETCOUNTERS] = 1;
1533 if (ir->nstlist == -1 && !bRerunMD)
1535 /* When bGStatEveryStep=FALSE, global_stat is only called
1536 * when we check the atom displacements, not at NS steps.
1537 * This means that also the bonded interaction count check is not
1538 * performed immediately after NS. Therefore a few MD steps could
1539 * be performed with missing interactions.
1540 * But wrong energies are never written to file,
1541 * since energies are only written after global_stat
1544 if (step >= nlh.step_nscheck)
1546 nlh.nabnsb = natoms_beyond_ns_buffer(ir,fr,&top->cgs,
1547 nlh.scale_tot,state->x);
1551 /* This is not necessarily true,
1552 * but step_nscheck is determined quite conservatively.
1558 /* In parallel we only have to check for checkpointing in steps
1559 * where we do global communication,
1560 * otherwise the other nodes don't know.
1562 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1565 run_time >= nchkpt*cpt_period*60.0)) &&
1566 gs.set[eglsCHKPT] == 0)
1568 gs.sig[eglsCHKPT] = 1;
1571 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1576 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1578 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1580 gmx_bool bIfRandomize;
1581 bIfRandomize = update_randomize_velocities(ir,step,mdatoms,state,upd,&top->idef,constr);
1582 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1583 if (constr && bIfRandomize)
1585 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1586 state,fr->bMolPBC,graph,f,
1587 &top->idef,tmp_vir,NULL,
1588 cr,nrnb,wcycle,upd,constr,
1589 bInitStep,TRUE,bCalcVir,vetanew);
1594 if (bIterativeCase && do_per_step(step,ir->nstpcouple))
1596 gmx_iterate_init(&iterate,TRUE);
1597 /* for iterations, we save these vectors, as we will be redoing the calculations */
1598 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1601 bFirstIterate = TRUE;
1602 while (bFirstIterate || iterate.bIterationActive)
1604 /* We now restore these vectors to redo the calculation with improved extended variables */
1605 if (iterate.bIterationActive)
1607 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1610 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1611 so scroll down for that logic */
1613 /* ######### START SECOND UPDATE STEP ################# */
1614 /* Box is changed in update() when we do pressure coupling,
1615 * but we should still use the old box for energy corrections and when
1616 * writing it to the energy file, so it matches the trajectory files for
1617 * the same timestep above. Make a copy in a separate array.
1619 copy_mat(state->box,lastbox);
1622 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1624 wallcycle_start(wcycle,ewcUPDATE);
1626 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1629 if (iterate.bIterationActive)
1637 /* we use a new value of scalevir to converge the iterations faster */
1638 scalevir = tracevir/trace(shake_vir);
1640 msmul(shake_vir,scalevir,shake_vir);
1641 m_add(force_vir,shake_vir,total_vir);
1642 clear_mat(shake_vir);
1644 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ3);
1645 /* We can only do Berendsen coupling after we have summed
1646 * the kinetic energy or virial. Since the happens
1647 * in global_state after update, we should only do it at
1648 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1653 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1654 update_pcouple(fplog,step,ir,state,pcoupl_mu,M,wcycle,
1660 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1662 /* velocity half-step update */
1663 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1664 bUpdateDoLR,fr->f_twin,fcd,
1665 ekind,M,wcycle,upd,FALSE,etrtVELOCITY2,
1666 cr,nrnb,constr,&top->idef);
1669 /* Above, initialize just copies ekinh into ekin,
1670 * it doesn't copy position (for VV),
1671 * and entire integrator for MD.
1676 copy_rvecn(state->x,cbuf,0,state->natoms);
1678 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1680 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1681 bUpdateDoLR,fr->f_twin,fcd,
1682 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1683 wallcycle_stop(wcycle,ewcUPDATE);
1685 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,state,
1686 fr->bMolPBC,graph,f,
1687 &top->idef,shake_vir,force_vir,
1688 cr,nrnb,wcycle,upd,constr,
1689 bInitStep,FALSE,bCalcVir,state->veta);
1693 /* erase F_EKIN and F_TEMP here? */
1694 /* just compute the kinetic energy at the half step to perform a trotter step */
1695 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1696 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1697 constr,NULL,FALSE,lastbox,
1698 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1699 cglo_flags | CGLO_TEMPERATURE
1701 wallcycle_start(wcycle,ewcUPDATE);
1702 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ4);
1703 /* now we know the scaling, we can compute the positions again again */
1704 copy_rvecn(cbuf,state->x,0,state->natoms);
1706 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1708 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1709 bUpdateDoLR,fr->f_twin,fcd,
1710 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1711 wallcycle_stop(wcycle,ewcUPDATE);
1713 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1714 /* are the small terms in the shake_vir here due
1715 * to numerical errors, or are they important
1716 * physically? I'm thinking they are just errors, but not completely sure.
1717 * For now, will call without actually constraining, constr=NULL*/
1718 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1719 state,fr->bMolPBC,graph,f,
1720 &top->idef,tmp_vir,force_vir,
1721 cr,nrnb,wcycle,upd,NULL,
1722 bInitStep,FALSE,bCalcVir,
1725 if (!bOK && !bFFscan)
1727 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1730 if (fr->bSepDVDL && fplog && do_log)
1732 fprintf(fplog,sepdvdlformat,"Constraint dV/dl",0.0,dvdl);
1734 enerd->term[F_DVDL_BONDED] += dvdl;
1738 /* Need to unshift here */
1739 unshift_self(graph,state->box,state->x);
1744 wallcycle_start(wcycle,ewcVSITECONSTR);
1747 shift_self(graph,state->box,state->x);
1749 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
1750 top->idef.iparams,top->idef.il,
1751 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
1755 unshift_self(graph,state->box,state->x);
1757 wallcycle_stop(wcycle,ewcVSITECONSTR);
1760 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1761 /* With Leap-Frog we can skip compute_globals at
1762 * non-communication steps, but we need to calculate
1763 * the kinetic energy one step before communication.
1765 if (bGStat || (!EI_VV(ir->eI)&&do_per_step(step+1,nstglobalcomm)))
1767 if (ir->nstlist == -1 && bFirstIterate)
1769 gs.sig[eglsNABNSB] = nlh.nabnsb;
1771 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1772 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1774 bFirstIterate ? &gs : NULL,
1775 (step_rel % gs.nstms == 0) &&
1776 (multisim_nsteps<0 || (step_rel<multisim_nsteps)),
1778 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1780 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1781 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1782 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1783 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1784 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1785 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1788 if (ir->nstlist == -1 && bFirstIterate)
1790 nlh.nabnsb = gs.set[eglsNABNSB];
1791 gs.set[eglsNABNSB] = 0;
1794 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1795 /* ############# END CALC EKIN AND PRESSURE ################# */
1797 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1798 the virial that should probably be addressed eventually. state->veta has better properies,
1799 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1800 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1802 if (iterate.bIterationActive &&
1803 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1804 trace(shake_vir),&tracevir))
1808 bFirstIterate = FALSE;
1811 /* only add constraint dvdl after constraints */
1812 enerd->term[F_DVDL_BONDED] += dvdl;
1813 if (!bVV || bRerunMD)
1815 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1816 sum_dhdl(enerd,state->lambda,ir->fepvals);
1818 update_box(fplog,step,ir,mdatoms,state,graph,f,
1819 ir->nstlist==-1 ? &nlh.scale_tot : NULL,pcoupl_mu,nrnb,wcycle,upd,bInitStep,FALSE);
1821 /* ################# END UPDATE STEP 2 ################# */
1822 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1824 /* The coordinates (x) were unshifted in update */
1825 if (bFFscan && (shellfc==NULL || bConverged))
1827 if (print_forcefield(fplog,enerd->term,mdatoms->homenr,
1829 &(top_global->mols),mdatoms->massT,pres))
1833 fprintf(stderr,"\n");
1839 /* We will not sum ekinh_old,
1840 * so signal that we still have to do it.
1842 bSumEkinhOld = TRUE;
1847 /* Only do GCT when the relaxation of shells (minimization) has converged,
1848 * otherwise we might be coupling to bogus energies.
1849 * In parallel we must always do this, because the other sims might
1853 /* Since this is called with the new coordinates state->x, I assume
1854 * we want the new box state->box too. / EL 20040121
1856 do_coupling(fplog,oenv,nfile,fnm,tcr,t,step,enerd->term,fr,
1858 mdatoms,&(top->idef),mu_aver,
1859 top_global->mols.nr,cr,
1860 state->box,total_vir,pres,
1861 mu_tot,state->x,f,bConverged);
1865 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1867 /* use the directly determined last velocity, not actually the averaged half steps */
1868 if (bTrotter && ir->eI==eiVV)
1870 enerd->term[F_EKIN] = last_ekin;
1872 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1876 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1880 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir,state,&MassQ);
1882 /* Check for excessively large energies */
1886 real etot_max = 1e200;
1888 real etot_max = 1e30;
1890 if (fabs(enerd->term[F_ETOT]) > etot_max)
1892 fprintf(stderr,"Energy too large (%g), giving up\n",
1893 enerd->term[F_ETOT]);
1896 /* ######### END PREPARING EDR OUTPUT ########### */
1898 /* Time for performance */
1899 if (((step % stepout) == 0) || bLastStep)
1901 runtime_upd_proc(runtime);
1907 gmx_bool do_dr,do_or;
1909 if (fplog && do_log && bDoExpanded)
1911 /* only needed if doing expanded ensemble */
1912 PrintFreeEnergyInfoToFile(fplog,ir->fepvals,ir->expandedvals,ir->bSimTemp?ir->simtempvals:NULL,
1913 &df_history,state->fep_state,ir->nstlog,step);
1915 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1919 upd_mdebin(mdebin,bDoDHDL, TRUE,
1920 t,mdatoms->tmass,enerd,state,
1921 ir->fepvals,ir->expandedvals,lastbox,
1922 shake_vir,force_vir,total_vir,pres,
1923 ekind,mu_tot,constr);
1927 upd_mdebin_step(mdebin);
1930 do_dr = do_per_step(step,ir->nstdisreout);
1931 do_or = do_per_step(step,ir->nstorireout);
1933 print_ebin(outf->fp_ene,do_ene,do_dr,do_or,do_log?fplog:NULL,
1935 eprNORMAL,bCompact,mdebin,fcd,groups,&(ir->opts));
1937 if (ir->ePull != epullNO)
1939 pull_print_output(ir->pull,step,t);
1942 if (do_per_step(step,ir->nstlog))
1944 if(fflush(fplog) != 0)
1946 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of disk space?");
1952 /* Have to do this part after outputting the logfile and the edr file */
1953 state->fep_state = lamnew;
1954 for (i=0;i<efptNR;i++)
1956 state_global->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
1959 /* Remaining runtime */
1960 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1964 fprintf(stderr,"\n");
1966 print_time(stderr,runtime,step,ir,cr);
1969 /* Replica exchange */
1971 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1972 do_per_step(step,repl_ex_nst))
1974 bExchanged = replica_exchange(fplog,cr,repl_ex,
1978 if (bExchanged && DOMAINDECOMP(cr))
1980 dd_partition_system(fplog,step,cr,TRUE,1,
1981 state_global,top_global,ir,
1982 state,&f,mdatoms,top,fr,
1983 vsite,shellfc,constr,
1990 bStartingFromCpt = FALSE;
1992 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1993 /* With all integrators, except VV, we need to retain the pressure
1994 * at the current step for coupling at the next step.
1996 if ((state->flags & (1<<estPRES_PREV)) &&
1998 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
2000 /* Store the pressure in t_state for pressure coupling
2001 * at the next MD step.
2003 copy_mat(pres,state->pres_prev);
2006 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2008 if ( (membed!=NULL) && (!bLastStep) )
2010 rescale_membed(step_rel,membed,state_global->x);
2017 /* read next frame from input trajectory */
2018 bNotLastFrame = read_next_frame(oenv,status,&rerun_fr);
2023 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
2027 if (!bRerunMD || !rerun_fr.bStep)
2029 /* increase the MD step number */
2034 cycles = wallcycle_stop(wcycle,ewcSTEP);
2035 if (DOMAINDECOMP(cr) && wcycle)
2037 dd_cycles_add(cr->dd,cycles,ddCyclStep);
2040 if (bPMETuneRunning || bPMETuneTry)
2042 /* PME grid + cut-off optimization with GPUs or PME nodes */
2044 /* Count the total cycles over the last steps */
2045 cycles_pmes += cycles;
2047 /* We can only switch cut-off at NS steps */
2048 if (step % ir->nstlist == 0)
2050 /* PME grid + cut-off optimization with GPUs or PME nodes */
2053 if (DDMASTER(cr->dd))
2055 /* PME node load is too high, start tuning */
2056 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2058 dd_bcast(cr->dd,sizeof(gmx_bool),&bPMETuneRunning);
2060 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2062 bPMETuneTry = FALSE;
2065 if (bPMETuneRunning)
2067 /* init_step might not be a multiple of nstlist,
2068 * but the first cycle is always skipped anyhow.
2071 pme_load_balance(pme_loadbal,cr,
2072 (bVerbose && MASTER(cr)) ? stderr : NULL,
2074 ir,state,cycles_pmes,
2075 fr->ic,fr->nbv,&fr->pmedata,
2078 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2079 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2080 fr->rlist = fr->ic->rlist;
2081 fr->rlistlong = fr->ic->rlistlong;
2082 fr->rcoulomb = fr->ic->rcoulomb;
2083 fr->rvdw = fr->ic->rvdw;
2089 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2090 gs.set[eglsRESETCOUNTERS] != 0)
2092 /* Reset all the counters related to performance over the run */
2093 reset_all_counters(fplog,cr,step,&step_rel,ir,wcycle,nrnb,runtime,
2094 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2095 wcycle_set_reset_counters(wcycle,-1);
2096 /* Correct max_hours for the elapsed time */
2097 max_hours -= run_time/(60.0*60.0);
2098 bResetCountersHalfMaxH = FALSE;
2099 gs.set[eglsRESETCOUNTERS] = 0;
2103 /* End of main MD loop */
2107 runtime_end(runtime);
2109 if (bRerunMD && MASTER(cr))
2114 if (!(cr->duty & DUTY_PME))
2116 /* Tell the PME only node to finish */
2117 gmx_pme_send_finish(cr);
2122 if (ir->nstcalcenergy > 0 && !bRerunMD)
2124 print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
2125 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
2133 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2135 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)));
2136 fprintf(fplog,"Average number of atoms that crossed the half buffer length: %.1f\n\n",nlh.ab/nlh.nns);
2139 if (pme_loadbal != NULL)
2141 pme_loadbal_done(pme_loadbal,fplog);
2144 if (shellfc && fplog)
2146 fprintf(fplog,"Fraction of iterations that converged: %.2f %%\n",
2147 (nconverged*100.0)/step_rel);
2148 fprintf(fplog,"Average number of force evaluations per MD step: %.2f\n\n",
2152 if (repl_ex_nst > 0 && MASTER(cr))
2154 print_replica_exchange_statistics(fplog,repl_ex);
2157 runtime->nsteps_done = step_rel;