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;
204 gmx_bool bUpdateDoLR;
207 atom_id *grpindex=NULL;
209 t_coupl_rec *tcr=NULL;
210 rvec *xcopy=NULL,*vcopy=NULL,*cbuf=NULL;
211 matrix boxcopy={{0}},lastbox;
213 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 bIterations = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
268 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
272 /* Since we don't know if the frames read are related in any way,
273 * rebuild the neighborlist at every step.
276 ir->nstcalcenergy = 1;
280 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
282 nstglobalcomm = check_nstglobalcomm(fplog,cr,nstglobalcomm,ir);
283 bGStatEveryStep = (nstglobalcomm == 1);
285 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
288 "To reduce the energy communication with nstlist = -1\n"
289 "the neighbor list validity should not be checked at every step,\n"
290 "this means that exact integration is not guaranteed.\n"
291 "The neighbor list validity is checked after:\n"
292 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
293 "In most cases this will result in exact integration.\n"
294 "This reduces the energy communication by a factor of 2 to 3.\n"
295 "If you want less energy communication, set nstlist > 3.\n\n");
298 if (bRerunMD || bFFscan)
302 groups = &top_global->groups;
305 init_md(fplog,cr,ir,oenv,&t,&t0,state_global->lambda,
306 &(state_global->fep_state),lam0,
307 nrnb,top_global,&upd,
308 nfile,fnm,&outf,&mdebin,
309 force_vir,shake_vir,mu_tot,&bSimAnn,&vcm,state_global,Flags);
311 clear_mat(total_vir);
313 /* Energy terms and groups */
315 init_enerdata(top_global->groups.grps[egcENER].nr,ir->fepvals->n_lambda,
317 if (DOMAINDECOMP(cr))
323 snew(f,top_global->natoms);
326 /* lambda Monte carlo random number generator */
329 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
331 /* copy the state into df_history */
332 copy_df_history(&df_history,&state_global->dfhist);
334 /* Kinetic energy data */
336 init_ekindata(fplog,top_global,&(ir->opts),ekind);
337 /* needed for iteration of constraints */
339 init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
340 /* Copy the cos acceleration to the groups struct */
341 ekind->cosacc.cos_accel = ir->cos_accel;
343 gstat = global_stat_init(ir);
346 /* Check for polarizable models and flexible constraints */
347 shellfc = init_shell_flexcon(fplog,
348 top_global,n_flexible_constraints(constr),
349 (ir->bContinuation ||
350 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
351 NULL : state_global->x);
355 #ifdef GMX_THREAD_MPI
356 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
358 set_deform_reference_box(upd,
359 deform_init_init_step_tpx,
360 deform_init_box_tpx);
361 #ifdef GMX_THREAD_MPI
362 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
367 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
368 if ((io > 2000) && MASTER(cr))
370 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
374 if (DOMAINDECOMP(cr)) {
375 top = dd_init_local_top(top_global);
378 dd_init_local_state(cr->dd,state_global,state);
380 if (DDMASTER(cr->dd) && ir->nstfout) {
381 snew(f_global,state_global->natoms);
385 /* Initialize the particle decomposition and split the topology */
386 top = split_system(fplog,top_global,ir,cr);
388 pd_cg_range(cr,&fr->cg0,&fr->hcg);
389 pd_at_range(cr,&a0,&a1);
391 top = gmx_mtop_generate_local_top(top_global,ir);
394 a1 = top_global->natoms;
397 forcerec_set_excl_load(fr,top,cr);
399 state = partdec_init_local_state(cr,state_global);
402 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
405 set_vsite_top(vsite,top,mdatoms,cr);
408 if (ir->ePBC != epbcNONE && !fr->bMolPBC) {
409 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
413 make_local_shells(cr,mdatoms,shellfc);
416 init_bonded_thread_force_reduction(fr,&top->idef);
418 if (ir->pull && PAR(cr)) {
419 dd_make_local_pull_groups(NULL,ir->pull,mdatoms);
423 if (DOMAINDECOMP(cr))
425 /* Distribute the charge groups over the nodes from the master node */
426 dd_partition_system(fplog,ir->init_step,cr,TRUE,1,
427 state_global,top_global,ir,
428 state,&f,mdatoms,top,fr,
429 vsite,shellfc,constr,
434 update_mdatoms(mdatoms,state->lambda[efptMASS]);
436 if (opt2bSet("-cpi",nfile,fnm))
438 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr);
442 bStateFromCP = FALSE;
449 /* Update mdebin with energy history if appending to output files */
450 if ( Flags & MD_APPENDFILES )
452 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
456 /* We might have read an energy history from checkpoint,
457 * free the allocated memory and reset the counts.
459 done_energyhistory(&state_global->enerhist);
460 init_energyhistory(&state_global->enerhist);
463 /* Set the initial energy history in state by updating once */
464 update_energyhistory(&state_global->enerhist,mdebin);
467 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
469 /* Set the random state if we read a checkpoint file */
470 set_stochd_state(upd,state);
473 if (state->flags & (1<<estMC_RNG))
475 set_mc_state(mcrng,state);
478 /* Initialize constraints */
480 if (!DOMAINDECOMP(cr))
481 set_constraints(constr,top,ir,mdatoms,cr);
484 /* Check whether we have to GCT stuff */
485 bTCR = ftp2bSet(efGCT,nfile,fnm);
488 fprintf(stderr,"Will do General Coupling Theory!\n");
490 gnx = top_global->mols.nr;
492 for(i=0; (i<gnx); i++) {
499 /* We need to be sure replica exchange can only occur
500 * when the energies are current */
501 check_nst_param(fplog,cr,"nstcalcenergy",ir->nstcalcenergy,
502 "repl_ex_nst",&repl_ex_nst);
503 /* This check needs to happen before inter-simulation
504 * signals are initialized, too */
506 if (repl_ex_nst > 0 && MASTER(cr))
508 repl_ex = init_replica_exchange(fplog,cr->ms,state_global,ir,
509 repl_ex_nst,repl_ex_nex,repl_ex_seed);
512 /* PME tuning is only supported with GPUs or PME nodes and not with rerun */
513 if ((Flags & MD_TUNEPME) &&
514 EEL_PME(fr->eeltype) &&
515 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
518 pme_loadbal_init(&pme_loadbal,ir,state->box,fr->ic,fr->pmedata);
520 if (cr->duty & DUTY_PME)
522 /* Start tuning right away, as we can't measure the load */
523 bPMETuneRunning = TRUE;
527 /* Separate PME nodes, we can measure the PP/PME load balance */
532 if (!ir->bContinuation && !bRerunMD)
534 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
536 /* Set the velocities of frozen particles to zero */
537 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
541 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
551 /* Constrain the initial coordinates and velocities */
552 do_constrain_first(fplog,constr,ir,mdatoms,state,f,
553 graph,cr,nrnb,fr,top,shake_vir);
557 /* Construct the virtual sites for the initial configuration */
558 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
559 top->idef.iparams,top->idef.il,
560 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
566 /* set free energy calculation frequency as the minimum of nstdhdl, nstexpanded, and nstrepl_ex_nst*/
567 nstfep = ir->fepvals->nstdhdl;
568 if (ir->bExpanded && (nstfep > ir->expandedvals->nstexpanded))
570 nstfep = ir->expandedvals->nstexpanded;
572 if (repl_ex_nst > 0 && repl_ex_nst > nstfep)
574 nstfep = repl_ex_nst;
577 /* I'm assuming we need global communication the first time! MRS */
578 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
579 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM:0)
580 | (bVV ? CGLO_PRESSURE:0)
581 | (bVV ? CGLO_CONSTRAINT:0)
582 | (bRerunMD ? CGLO_RERUNMD:0)
583 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN:0));
585 bSumEkinhOld = FALSE;
586 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
587 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
588 constr,NULL,FALSE,state->box,
589 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,cglo_flags);
590 if (ir->eI == eiVVAK) {
591 /* a second call to get the half step temperature initialized as well */
592 /* we do the same call as above, but turn the pressure off -- internally to
593 compute_globals, this is recognized as a velocity verlet half-step
594 kinetic energy calculation. This minimized excess variables, but
595 perhaps loses some logic?*/
597 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
598 NULL,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
599 constr,NULL,FALSE,state->box,
600 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
601 cglo_flags &~ (CGLO_STOPCM | CGLO_PRESSURE));
604 /* Calculate the initial half step temperature, and save the ekinh_old */
605 if (!(Flags & MD_STARTFROMCPT))
607 for(i=0; (i<ir->opts.ngtc); i++)
609 copy_mat(ekind->tcstat[i].ekinh,ekind->tcstat[i].ekinh_old);
614 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
615 and there is no previous step */
618 /* if using an iterative algorithm, we need to create a working directory for the state. */
621 bufstate = init_bufstate(state);
625 snew(xcopy,state->natoms);
626 snew(vcopy,state->natoms);
627 copy_rvecn(state->x,xcopy,0,state->natoms);
628 copy_rvecn(state->v,vcopy,0,state->natoms);
629 copy_mat(state->box,boxcopy);
632 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
633 temperature control */
634 trotter_seq = init_npt_vars(ir,state,&MassQ,bTrotter);
638 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
641 "RMS relative constraint deviation after constraining: %.2e\n",
642 constr_rmsd(constr,FALSE));
644 if (EI_STATE_VELOCITY(ir->eI))
646 fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
650 fprintf(stderr,"starting md rerun '%s', reading coordinates from"
651 " input trajectory '%s'\n\n",
652 *(top_global->name),opt2fn("-rerun",nfile,fnm));
655 fprintf(stderr,"Calculated time to finish depends on nsteps from "
656 "run input file,\nwhich may not correspond to the time "
657 "needed to process input trajectory.\n\n");
663 fprintf(stderr,"starting mdrun '%s'\n",
664 *(top_global->name));
667 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
671 sprintf(tbuf,"%s","infinite");
673 if (ir->init_step > 0)
675 fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
676 gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
677 gmx_step_str(ir->init_step,sbuf2),
678 ir->init_step*ir->delta_t);
682 fprintf(stderr,"%s steps, %s ps.\n",
683 gmx_step_str(ir->nsteps,sbuf),tbuf);
689 /* Set and write start time */
690 runtime_start(runtime);
691 print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
692 wallcycle_start(wcycle,ewcRUN);
698 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
700 chkpt_ret=fcCheckPointParallel( cr->nodeid,
702 if ( chkpt_ret == 0 )
703 gmx_fatal( 3,__FILE__,__LINE__, "Checkpoint error on step %d\n", 0 );
707 /***********************************************************
711 ************************************************************/
713 /* if rerunMD then read coordinates and velocities from input trajectory */
716 if (getenv("GMX_FORCE_UPDATE"))
724 bNotLastFrame = read_first_frame(oenv,&status,
725 opt2fn("-rerun",nfile,fnm),
726 &rerun_fr,TRX_NEED_X | TRX_READ_V);
727 if (rerun_fr.natoms != top_global->natoms)
730 "Number of atoms in trajectory (%d) does not match the "
731 "run input file (%d)\n",
732 rerun_fr.natoms,top_global->natoms);
734 if (ir->ePBC != epbcNONE)
738 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);
740 if (max_cutoff2(ir->ePBC,rerun_fr.box) < sqr(fr->rlistlong))
742 gmx_fatal(FARGS,"Rerun trajectory frame step %d time %f has too small box dimensions",rerun_fr.step,rerun_fr.time);
749 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
752 if (ir->ePBC != epbcNONE)
754 /* Set the shift vectors.
755 * Necessary here when have a static box different from the tpr box.
757 calc_shifts(rerun_fr.box,fr->shift_vec);
761 /* loop over MD steps or if rerunMD to end of input trajectory */
763 /* Skip the first Nose-Hoover integration when we get the state from tpx */
764 bStateFromTPX = !bStateFromCP;
765 bInitStep = bFirstStep && (bStateFromTPX || bVV);
766 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
768 bSumEkinhOld = FALSE;
771 init_global_signals(&gs,cr,ir,repl_ex_nst);
773 step = ir->init_step;
776 if (ir->nstlist == -1)
778 init_nlistheuristics(&nlh,bGStatEveryStep,step);
781 if (MULTISIM(cr) && (repl_ex_nst <=0 ))
783 /* check how many steps are left in other sims */
784 multisim_nsteps=get_multisim_nsteps(cr, ir->nsteps);
788 /* and stop now if we should */
789 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
790 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
791 while (!bLastStep || (bRerunMD && bNotLastFrame)) {
793 wallcycle_start(wcycle,ewcSTEP);
796 if (rerun_fr.bStep) {
797 step = rerun_fr.step;
798 step_rel = step - ir->init_step;
800 if (rerun_fr.bTime) {
810 bLastStep = (step_rel == ir->nsteps);
811 t = t0 + step*ir->delta_t;
814 if (ir->efep != efepNO || ir->bSimTemp)
816 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
817 requiring different logic. */
819 set_current_lambdas(step,ir->fepvals,bRerunMD,&rerun_fr,state_global,state,lam0);
820 bDoDHDL = do_per_step(step,ir->fepvals->nstdhdl);
821 bDoFEP = (do_per_step(step,nstfep) && (ir->efep != efepNO));
822 bDoExpanded = (do_per_step(step,ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
827 update_annealing_target_temp(&(ir->opts),t);
832 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
834 for(i=0; i<state_global->natoms; i++)
836 copy_rvec(rerun_fr.x[i],state_global->x[i]);
840 for(i=0; i<state_global->natoms; i++)
842 copy_rvec(rerun_fr.v[i],state_global->v[i]);
847 for(i=0; i<state_global->natoms; i++)
849 clear_rvec(state_global->v[i]);
853 fprintf(stderr,"\nWARNING: Some frames do not contain velocities.\n"
854 " Ekin, temperature and pressure are incorrect,\n"
855 " the virial will be incorrect when constraints are present.\n"
857 bRerunWarnNoV = FALSE;
861 copy_mat(rerun_fr.box,state_global->box);
862 copy_mat(state_global->box,state->box);
864 if (vsite && (Flags & MD_RERUN_VSITE))
866 if (DOMAINDECOMP(cr))
868 gmx_fatal(FARGS,"Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
872 /* Following is necessary because the graph may get out of sync
873 * with the coordinates if we only have every N'th coordinate set
875 mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
876 shift_self(graph,state->box,state->x);
878 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
879 top->idef.iparams,top->idef.il,
880 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
883 unshift_self(graph,state->box,state->x);
888 /* Stop Center of Mass motion */
889 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step,ir->nstcomm));
891 /* Copy back starting coordinates in case we're doing a forcefield scan */
894 for(ii=0; (ii<state->natoms); ii++)
896 copy_rvec(xcopy[ii],state->x[ii]);
897 copy_rvec(vcopy[ii],state->v[ii]);
899 copy_mat(boxcopy,state->box);
904 /* for rerun MD always do Neighbour Searching */
905 bNS = (bFirstStep || ir->nstlist != 0);
910 /* Determine whether or not to do Neighbour Searching and LR */
911 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
913 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
914 (ir->nstlist == -1 && nlh.nabnsb > 0));
916 if (bNS && ir->nstlist == -1)
918 set_nlistheuristics(&nlh,bFirstStep || bExchanged || bDoFEP, step);
922 /* check whether we should stop because another simulation has
926 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
927 (multisim_nsteps != ir->nsteps) )
934 "Stopping simulation %d because another one has finished\n",
938 gs.sig[eglsCHKPT] = 1;
943 /* < 0 means stop at next step, > 0 means stop at next NS step */
944 if ( (gs.set[eglsSTOPCOND] < 0 ) ||
945 ( (gs.set[eglsSTOPCOND] > 0 ) && ( bNS || ir->nstlist==0)) )
950 /* Determine whether or not to update the Born radii if doing GB */
951 bBornRadii=bFirstStep;
952 if (ir->implicit_solvent && (step % ir->nstgbradii==0))
957 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
958 do_verbose = bVerbose &&
959 (step % stepout == 0 || bFirstStep || bLastStep);
961 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
969 bMasterState = FALSE;
970 /* Correct the new box if it is too skewed */
971 if (DYNAMIC_BOX(*ir))
973 if (correct_box(fplog,step,state->box,graph))
978 if (DOMAINDECOMP(cr) && bMasterState)
980 dd_collect_state(cr->dd,state,state_global);
984 if (DOMAINDECOMP(cr))
986 /* Repartition the domain decomposition */
987 wallcycle_start(wcycle,ewcDOMDEC);
988 dd_partition_system(fplog,step,cr,
989 bMasterState,nstglobalcomm,
990 state_global,top_global,ir,
991 state,&f,mdatoms,top,fr,
992 vsite,shellfc,constr,
994 do_verbose && !bPMETuneRunning);
995 wallcycle_stop(wcycle,ewcDOMDEC);
996 /* If using an iterative integrator, reallocate space to match the decomposition */
1000 if (MASTER(cr) && do_log && !bFFscan)
1002 print_ebin_header(fplog,step,t,state->lambda[efptFEP]); /* can we improve the information printed here? */
1005 if (ir->efep != efepNO)
1007 update_mdatoms(mdatoms,state->lambda[efptMASS]);
1010 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1013 /* We need the kinetic energy at minus the half step for determining
1014 * the full step kinetic energy and possibly for T-coupling.*/
1015 /* This may not be quite working correctly yet . . . . */
1016 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1017 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1018 constr,NULL,FALSE,state->box,
1019 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1020 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1022 clear_mat(force_vir);
1024 /* Ionize the atoms if necessary */
1027 ionize(fplog,oenv,mdatoms,top_global,t,ir,state->x,state->v,
1028 mdatoms->start,mdatoms->start+mdatoms->homenr,state->box,cr);
1031 /* Update force field in ffscan program */
1034 if (update_forcefield(fplog,
1036 mdatoms->nr,state->x,state->box))
1044 /* We write a checkpoint at this MD step when:
1045 * either at an NS step when we signalled through gs,
1046 * or at the last step (but not when we do not want confout),
1047 * but never at the first step or with rerun.
1049 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1050 (bLastStep && (Flags & MD_CONFOUT))) &&
1051 step > ir->init_step && !bRerunMD);
1054 gs.set[eglsCHKPT] = 0;
1057 /* Determine the energy and pressure:
1058 * at nstcalcenergy steps and at energy output steps (set below).
1060 if (EI_VV(ir->eI) && (!bInitStep))
1062 /* for vv, the first half actually corresponds to the last step */
1063 bCalcEner = do_per_step(step-1,ir->nstcalcenergy);
1067 bCalcEner = do_per_step(step,ir->nstcalcenergy);
1069 bCalcVir = bCalcEner ||
1070 (ir->epc != epcNO && do_per_step(step,ir->nstpcouple));
1072 /* Do we need global communication ? */
1073 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1074 do_per_step(step,nstglobalcomm) ||
1075 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1077 do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
1079 if (do_ene || do_log)
1086 /* these CGLO_ options remain the same throughout the iteration */
1087 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1088 (bGStat ? CGLO_GSTAT : 0)
1091 force_flags = (GMX_FORCE_STATECHANGED |
1092 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1093 GMX_FORCE_ALLFORCES |
1095 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1096 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1097 (bDoFEP ? GMX_FORCE_DHDL : 0)
1102 if(do_per_step(step,ir->nstcalclr))
1104 force_flags |= GMX_FORCE_DO_LR;
1110 /* Now is the time to relax the shells */
1111 count=relax_shell_flexcon(fplog,cr,bVerbose,bFFscan ? step+1 : step,
1113 bStopCM,top,top_global,
1115 state,f,force_vir,mdatoms,
1116 nrnb,wcycle,graph,groups,
1117 shellfc,fr,bBornRadii,t,mu_tot,
1118 state->natoms,&bConverged,vsite,
1129 /* The coordinates (x) are shifted (to get whole molecules)
1131 * This is parallellized as well, and does communication too.
1132 * Check comments in sim_util.c
1134 do_force(fplog,cr,ir,step,nrnb,wcycle,top,top_global,groups,
1135 state->box,state->x,&state->hist,
1136 f,force_vir,mdatoms,enerd,fcd,
1137 state->lambda,graph,
1138 fr,vsite,mu_tot,t,outf->fp_field,ed,bBornRadii,
1139 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1144 mu_aver = calc_mu_aver(cr,state->x,mdatoms->chargeA,
1145 mu_tot,&top_global->mols,mdatoms,gnx,grpindex);
1148 if (bTCR && bFirstStep)
1150 tcr=init_coupling(fplog,nfile,fnm,cr,fr,mdatoms,&(top->idef));
1151 fprintf(fplog,"Done init_coupling\n");
1155 if (bVV && !bStartingFromCpt && !bRerunMD)
1156 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1158 if (ir->eI==eiVV && bInitStep)
1160 /* if using velocity verlet with full time step Ekin,
1161 * take the first half step only to compute the
1162 * virial for the first step. From there,
1163 * revert back to the initial coordinates
1164 * so that the input is actually the initial step.
1166 copy_rvecn(state->v,cbuf,0,state->natoms); /* should make this better for parallelizing? */
1168 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1169 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ1);
1172 /* If we are using twin-range interactions where the long-range component
1173 * is only evaluated every nstcalclr>1 steps, we should do a special update
1174 * step to combine the long-range forces on these steps.
1175 * For nstcalclr=1 this is not done, since the forces would have been added
1176 * directly to the short-range forces already.
1178 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1180 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,
1181 f,bUpdateDoLR,fr->f_twin,fcd,
1182 ekind,M,wcycle,upd,bInitStep,etrtVELOCITY1,
1183 cr,nrnb,constr,&top->idef);
1187 gmx_iterate_init(&iterate,bIterations && !bInitStep);
1189 /* for iterations, we save these vectors, as we will be self-consistently iterating
1192 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1194 /* save the state */
1195 if (bIterations && iterate.bIterate) {
1196 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1199 bFirstIterate = TRUE;
1200 while (bFirstIterate || (bIterations && iterate.bIterate))
1202 if (bIterations && iterate.bIterate)
1204 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1205 if (bFirstIterate && bTrotter)
1207 /* The first time through, we need a decent first estimate
1208 of veta(t+dt) to compute the constraints. Do
1209 this by computing the box volume part of the
1210 trotter integration at this time. Nothing else
1211 should be changed by this routine here. If
1212 !(first time), we start with the previous value
1215 veta_save = state->veta;
1216 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ0);
1217 vetanew = state->veta;
1218 state->veta = veta_save;
1223 if ( !bRerunMD || rerun_fr.bV || bForceUpdate) { /* Why is rerun_fr.bV here? Unclear. */
1226 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1227 state,fr->bMolPBC,graph,f,
1228 &top->idef,shake_vir,NULL,
1229 cr,nrnb,wcycle,upd,constr,
1230 bInitStep,TRUE,bCalcVir,vetanew);
1232 if (!bOK && !bFFscan)
1234 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1239 { /* Need to unshift here if a do_force has been
1240 called in the previous step */
1241 unshift_self(graph,state->box,state->x);
1245 /* if VV, compute the pressure and constraints */
1246 /* For VV2, we strictly only need this if using pressure
1247 * control, but we really would like to have accurate pressures
1249 * Think about ways around this in the future?
1250 * For now, keep this choice in comments.
1252 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1253 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1255 bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK));
1256 if (bCalcEner && ir->eI==eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1258 bSumEkinhOld = TRUE;
1260 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1261 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1262 constr,NULL,FALSE,state->box,
1263 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1266 | (bStopCM ? CGLO_STOPCM : 0)
1267 | (bTemp ? CGLO_TEMPERATURE:0)
1268 | (bPres ? CGLO_PRESSURE : 0)
1269 | (bPres ? CGLO_CONSTRAINT : 0)
1270 | ((bIterations && iterate.bIterate) ? CGLO_ITERATE : 0)
1271 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1274 /* explanation of above:
1275 a) We compute Ekin at the full time step
1276 if 1) we are using the AveVel Ekin, and it's not the
1277 initial step, or 2) if we are using AveEkin, but need the full
1278 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1279 b) If we are using EkinAveEkin for the kinetic energy for the temperture control, we still feed in
1280 EkinAveVel because it's needed for the pressure */
1282 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1287 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ2);
1294 /* We need the kinetic energy at minus the half step for determining
1295 * the full step kinetic energy and possibly for T-coupling.*/
1296 /* This may not be quite working correctly yet . . . . */
1297 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1298 wcycle,enerd,NULL,NULL,NULL,NULL,mu_tot,
1299 constr,NULL,FALSE,state->box,
1300 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1301 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1305 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1310 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1311 state->veta,&vetanew))
1315 bFirstIterate = FALSE;
1318 if (bTrotter && !bInitStep) {
1319 enerd->term[F_DVDL_BONDED] += dvdl; /* only add after iterations */
1320 copy_mat(shake_vir,state->svir_prev);
1321 copy_mat(force_vir,state->fvir_prev);
1322 if (IR_NVT_TROTTER(ir) && ir->eI==eiVV) {
1323 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1324 enerd->term[F_TEMP] = sum_ekin(&(ir->opts),ekind,NULL,(ir->eI==eiVV),FALSE,FALSE);
1325 enerd->term[F_EKIN] = trace(ekind->ekin);
1328 /* if it's the initial step, we performed this first step just to get the constraint virial */
1329 if (bInitStep && ir->eI==eiVV) {
1330 copy_rvecn(cbuf,state->v,0,state->natoms);
1333 if (fr->bSepDVDL && fplog && do_log)
1335 fprintf(fplog,sepdvdlformat,"Constraint",0.0,dvdl);
1337 enerd->term[F_DVDL_BONDED] += dvdl;
1340 /* MRS -- now done iterating -- compute the conserved quantity */
1342 saved_conserved_quantity = compute_conserved_from_auxiliary(ir,state,&MassQ);
1345 last_ekin = enerd->term[F_EKIN];
1347 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1349 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1351 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1352 sum_dhdl(enerd,state->lambda,ir->fepvals);
1355 /* ######## END FIRST UPDATE STEP ############## */
1356 /* ######## If doing VV, we now have v(dt) ###### */
1358 /* perform extended ensemble sampling in lambda - we don't
1359 actually move to the new state before outputting
1360 statistics, but if performing simulated tempering, we
1361 do update the velocities and the tau_t. */
1363 lamnew = ExpandedEnsembleDynamics(fplog,ir,enerd,state,&MassQ,&df_history,step,mcrng,state->v,mdatoms);
1365 /* ################## START TRAJECTORY OUTPUT ################# */
1367 /* Now we have the energies and forces corresponding to the
1368 * coordinates at time t. We must output all of this before
1370 * for RerunMD t is read from input trajectory
1373 if (do_per_step(step,ir->nstxout)) { mdof_flags |= MDOF_X; }
1374 if (do_per_step(step,ir->nstvout)) { mdof_flags |= MDOF_V; }
1375 if (do_per_step(step,ir->nstfout)) { mdof_flags |= MDOF_F; }
1376 if (do_per_step(step,ir->nstxtcout)) { mdof_flags |= MDOF_XTC; }
1377 if (bCPT) { mdof_flags |= MDOF_CPT; };
1379 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1382 /* Enforce writing positions and velocities at end of run */
1383 mdof_flags |= (MDOF_X | MDOF_V);
1388 fcReportProgress( ir->nsteps, step );
1390 /* sync bCPT and fc record-keeping */
1391 if (bCPT && MASTER(cr))
1392 fcRequestCheckPoint();
1395 if (mdof_flags != 0)
1397 wallcycle_start(wcycle,ewcTRAJ);
1400 if (state->flags & (1<<estLD_RNG))
1402 get_stochd_state(upd,state);
1404 if (state->flags & (1<<estMC_RNG))
1406 get_mc_state(mcrng,state);
1412 state_global->ekinstate.bUpToDate = FALSE;
1416 update_ekinstate(&state_global->ekinstate,ekind);
1417 state_global->ekinstate.bUpToDate = TRUE;
1419 update_energyhistory(&state_global->enerhist,mdebin);
1420 if (ir->efep!=efepNO || ir->bSimTemp)
1422 state_global->fep_state = state->fep_state; /* MRS: seems kludgy. The code should be
1423 structured so this isn't necessary.
1424 Note this reassignment is only necessary
1425 for single threads.*/
1426 copy_df_history(&state_global->dfhist,&df_history);
1430 write_traj(fplog,cr,outf,mdof_flags,top_global,
1431 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
1438 if (bLastStep && step_rel == ir->nsteps &&
1439 (Flags & MD_CONFOUT) && MASTER(cr) &&
1440 !bRerunMD && !bFFscan)
1442 /* x and v have been collected in write_traj,
1443 * because a checkpoint file will always be written
1446 fprintf(stderr,"\nWriting final coordinates.\n");
1449 /* Make molecules whole only for confout writing */
1450 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
1452 write_sto_conf_mtop(ftp2fn(efSTO,nfile,fnm),
1453 *top_global->name,top_global,
1454 state_global->x,state_global->v,
1455 ir->ePBC,state->box);
1458 wallcycle_stop(wcycle,ewcTRAJ);
1461 /* kludge -- virial is lost with restart for NPT control. Must restart */
1462 if (bStartingFromCpt && bVV)
1464 copy_mat(state->svir_prev,shake_vir);
1465 copy_mat(state->fvir_prev,force_vir);
1467 /* ################## END TRAJECTORY OUTPUT ################ */
1469 /* Determine the wallclock run time up till now */
1470 run_time = gmx_gettime() - (double)runtime->real;
1472 /* Check whether everything is still allright */
1473 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1474 #ifdef GMX_THREAD_MPI
1479 /* this is just make gs.sig compatible with the hack
1480 of sending signals around by MPI_Reduce with together with
1482 if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
1483 gs.sig[eglsSTOPCOND]=1;
1484 if ( gmx_get_stop_condition() == gmx_stop_cond_next )
1485 gs.sig[eglsSTOPCOND]=-1;
1486 /* < 0 means stop at next step, > 0 means stop at next NS step */
1490 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1491 gmx_get_signal_name(),
1492 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1496 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1497 gmx_get_signal_name(),
1498 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
1500 handled_stop_condition=(int)gmx_get_stop_condition();
1502 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1503 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1504 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1506 /* Signal to terminate the run */
1507 gs.sig[eglsSTOPCOND] = 1;
1510 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1512 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
1515 if (bResetCountersHalfMaxH && MASTER(cr) &&
1516 run_time > max_hours*60.0*60.0*0.495)
1518 gs.sig[eglsRESETCOUNTERS] = 1;
1521 if (ir->nstlist == -1 && !bRerunMD)
1523 /* When bGStatEveryStep=FALSE, global_stat is only called
1524 * when we check the atom displacements, not at NS steps.
1525 * This means that also the bonded interaction count check is not
1526 * performed immediately after NS. Therefore a few MD steps could
1527 * be performed with missing interactions.
1528 * But wrong energies are never written to file,
1529 * since energies are only written after global_stat
1532 if (step >= nlh.step_nscheck)
1534 nlh.nabnsb = natoms_beyond_ns_buffer(ir,fr,&top->cgs,
1535 nlh.scale_tot,state->x);
1539 /* This is not necessarily true,
1540 * but step_nscheck is determined quite conservatively.
1546 /* In parallel we only have to check for checkpointing in steps
1547 * where we do global communication,
1548 * otherwise the other nodes don't know.
1550 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1553 run_time >= nchkpt*cpt_period*60.0)) &&
1554 gs.set[eglsCHKPT] == 0)
1556 gs.sig[eglsCHKPT] = 1;
1560 /* at the start of step, randomize the velocities */
1561 if (ETC_ANDERSEN(ir->etc) && EI_VV(ir->eI))
1563 gmx_bool bDoAndersenConstr;
1564 bDoAndersenConstr = (constr && update_randomize_velocities(ir,step,mdatoms,state,upd,&top->idef,constr));
1565 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1566 if (bDoAndersenConstr)
1568 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1569 state,fr->bMolPBC,graph,f,
1570 &top->idef,tmp_vir,NULL,
1571 cr,nrnb,wcycle,upd,constr,
1572 bInitStep,TRUE,bCalcVir,vetanew);
1578 gmx_iterate_init(&iterate,bIterations);
1581 /* for iterations, we save these vectors, as we will be redoing the calculations */
1582 if (bIterations && iterate.bIterate)
1584 copy_coupling_state(state,bufstate,ekind,ekind_save,&(ir->opts));
1586 bFirstIterate = TRUE;
1587 while (bFirstIterate || (bIterations && iterate.bIterate))
1589 /* We now restore these vectors to redo the calculation with improved extended variables */
1592 copy_coupling_state(bufstate,state,ekind_save,ekind,&(ir->opts));
1595 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1596 so scroll down for that logic */
1598 /* ######### START SECOND UPDATE STEP ################# */
1599 /* Box is changed in update() when we do pressure coupling,
1600 * but we should still use the old box for energy corrections and when
1601 * writing it to the energy file, so it matches the trajectory files for
1602 * the same timestep above. Make a copy in a separate array.
1604 copy_mat(state->box,lastbox);
1607 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1609 wallcycle_start(wcycle,ewcUPDATE);
1611 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1614 if (bIterations && iterate.bIterate)
1622 /* we use a new value of scalevir to converge the iterations faster */
1623 scalevir = tracevir/trace(shake_vir);
1625 msmul(shake_vir,scalevir,shake_vir);
1626 m_add(force_vir,shake_vir,total_vir);
1627 clear_mat(shake_vir);
1629 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ3);
1630 /* We can only do Berendsen coupling after we have summed
1631 * the kinetic energy or virial. Since the happens
1632 * in global_state after update, we should only do it at
1633 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1638 update_tcouple(fplog,step,ir,state,ekind,wcycle,upd,&MassQ,mdatoms);
1639 update_pcouple(fplog,step,ir,state,pcoupl_mu,M,wcycle,
1645 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1647 /* velocity half-step update */
1648 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1649 bUpdateDoLR,fr->f_twin,fcd,
1650 ekind,M,wcycle,upd,FALSE,etrtVELOCITY2,
1651 cr,nrnb,constr,&top->idef);
1654 /* Above, initialize just copies ekinh into ekin,
1655 * it doesn't copy position (for VV),
1656 * and entire integrator for MD.
1661 copy_rvecn(state->x,cbuf,0,state->natoms);
1663 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1665 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1666 bUpdateDoLR,fr->f_twin,fcd,
1667 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1668 wallcycle_stop(wcycle,ewcUPDATE);
1670 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,state,
1671 fr->bMolPBC,graph,f,
1672 &top->idef,shake_vir,force_vir,
1673 cr,nrnb,wcycle,upd,constr,
1674 bInitStep,FALSE,bCalcVir,state->veta);
1678 /* erase F_EKIN and F_TEMP here? */
1679 /* just compute the kinetic energy at the half step to perform a trotter step */
1680 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1681 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1682 constr,NULL,FALSE,lastbox,
1683 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1684 cglo_flags | CGLO_TEMPERATURE
1686 wallcycle_start(wcycle,ewcUPDATE);
1687 trotter_update(ir,step,ekind,enerd,state,total_vir,mdatoms,&MassQ,trotter_seq,ettTSEQ4);
1688 /* now we know the scaling, we can compute the positions again again */
1689 copy_rvecn(cbuf,state->x,0,state->natoms);
1691 bUpdateDoLR = (fr->bTwinRange && do_per_step(step,ir->nstcalclr));
1693 update_coords(fplog,step,ir,mdatoms,state,fr->bMolPBC,f,
1694 bUpdateDoLR,fr->f_twin,fcd,
1695 ekind,M,wcycle,upd,bInitStep,etrtPOSITION,cr,nrnb,constr,&top->idef);
1696 wallcycle_stop(wcycle,ewcUPDATE);
1698 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1699 /* are the small terms in the shake_vir here due
1700 * to numerical errors, or are they important
1701 * physically? I'm thinking they are just errors, but not completely sure.
1702 * For now, will call without actually constraining, constr=NULL*/
1703 update_constraints(fplog,step,&dvdl,ir,ekind,mdatoms,
1704 state,fr->bMolPBC,graph,f,
1705 &top->idef,tmp_vir,force_vir,
1706 cr,nrnb,wcycle,upd,NULL,
1707 bInitStep,FALSE,bCalcVir,
1710 if (!bOK && !bFFscan)
1712 gmx_fatal(FARGS,"Constraint error: Shake, Lincs or Settle could not solve the constrains");
1715 if (fr->bSepDVDL && fplog && do_log)
1717 fprintf(fplog,sepdvdlformat,"Constraint dV/dl",0.0,dvdl);
1719 enerd->term[F_DVDL_BONDED] += dvdl;
1723 /* Need to unshift here */
1724 unshift_self(graph,state->box,state->x);
1729 wallcycle_start(wcycle,ewcVSITECONSTR);
1732 shift_self(graph,state->box,state->x);
1734 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,state->v,
1735 top->idef.iparams,top->idef.il,
1736 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
1740 unshift_self(graph,state->box,state->x);
1742 wallcycle_stop(wcycle,ewcVSITECONSTR);
1745 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1746 /* With Leap-Frog we can skip compute_globals at
1747 * non-communication steps, but we need to calculate
1748 * the kinetic energy one step before communication.
1750 if (bGStat || do_per_step(step+1,nstglobalcomm) ||
1753 if (ir->nstlist == -1 && bFirstIterate)
1755 gs.sig[eglsNABNSB] = nlh.nabnsb;
1757 compute_globals(fplog,gstat,cr,ir,fr,ekind,state,state_global,mdatoms,nrnb,vcm,
1758 wcycle,enerd,force_vir,shake_vir,total_vir,pres,mu_tot,
1760 bFirstIterate ? &gs : NULL,
1761 (step_rel % gs.nstms == 0) &&
1762 (multisim_nsteps<0 || (step_rel<multisim_nsteps)),
1764 top_global,&pcurr,top_global->natoms,&bSumEkinhOld,
1766 | (!EI_VV(ir->eI) ? CGLO_ENERGY : 0)
1767 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1768 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1769 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1770 | (bIterations && iterate.bIterate ? CGLO_ITERATE : 0)
1771 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1774 if (ir->nstlist == -1 && bFirstIterate)
1776 nlh.nabnsb = gs.set[eglsNABNSB];
1777 gs.set[eglsNABNSB] = 0;
1780 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1781 /* ############# END CALC EKIN AND PRESSURE ################# */
1783 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1784 the virial that should probably be addressed eventually. state->veta has better properies,
1785 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1786 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1789 done_iterating(cr,fplog,step,&iterate,bFirstIterate,
1790 trace(shake_vir),&tracevir))
1794 bFirstIterate = FALSE;
1797 /* only add constraint dvdl after constraints */
1798 enerd->term[F_DVDL_BONDED] += dvdl;
1801 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1802 sum_dhdl(enerd,state->lambda,ir->fepvals);
1804 update_box(fplog,step,ir,mdatoms,state,graph,f,
1805 ir->nstlist==-1 ? &nlh.scale_tot : NULL,pcoupl_mu,nrnb,wcycle,upd,bInitStep,FALSE);
1807 /* ################# END UPDATE STEP 2 ################# */
1808 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1810 /* The coordinates (x) were unshifted in update */
1811 if (bFFscan && (shellfc==NULL || bConverged))
1813 if (print_forcefield(fplog,enerd->term,mdatoms->homenr,
1815 &(top_global->mols),mdatoms->massT,pres))
1819 fprintf(stderr,"\n");
1825 /* We will not sum ekinh_old,
1826 * so signal that we still have to do it.
1828 bSumEkinhOld = TRUE;
1833 /* Only do GCT when the relaxation of shells (minimization) has converged,
1834 * otherwise we might be coupling to bogus energies.
1835 * In parallel we must always do this, because the other sims might
1839 /* Since this is called with the new coordinates state->x, I assume
1840 * we want the new box state->box too. / EL 20040121
1842 do_coupling(fplog,oenv,nfile,fnm,tcr,t,step,enerd->term,fr,
1844 mdatoms,&(top->idef),mu_aver,
1845 top_global->mols.nr,cr,
1846 state->box,total_vir,pres,
1847 mu_tot,state->x,f,bConverged);
1851 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1853 /* use the directly determined last velocity, not actually the averaged half steps */
1854 if (bTrotter && ir->eI==eiVV)
1856 enerd->term[F_EKIN] = last_ekin;
1858 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1862 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1866 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir,state,&MassQ);
1868 /* Check for excessively large energies */
1872 real etot_max = 1e200;
1874 real etot_max = 1e30;
1876 if (fabs(enerd->term[F_ETOT]) > etot_max)
1878 fprintf(stderr,"Energy too large (%g), giving up\n",
1879 enerd->term[F_ETOT]);
1882 /* ######### END PREPARING EDR OUTPUT ########### */
1884 /* Time for performance */
1885 if (((step % stepout) == 0) || bLastStep)
1887 runtime_upd_proc(runtime);
1893 gmx_bool do_dr,do_or;
1895 if (fplog && do_log && bDoExpanded)
1897 /* only needed if doing expanded ensemble */
1898 PrintFreeEnergyInfoToFile(fplog,ir->fepvals,ir->expandedvals,ir->bSimTemp?ir->simtempvals:NULL,
1899 &df_history,state->fep_state,ir->nstlog,step);
1901 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1905 upd_mdebin(mdebin,bDoDHDL, TRUE,
1906 t,mdatoms->tmass,enerd,state,
1907 ir->fepvals,ir->expandedvals,lastbox,
1908 shake_vir,force_vir,total_vir,pres,
1909 ekind,mu_tot,constr);
1913 upd_mdebin_step(mdebin);
1916 do_dr = do_per_step(step,ir->nstdisreout);
1917 do_or = do_per_step(step,ir->nstorireout);
1919 print_ebin(outf->fp_ene,do_ene,do_dr,do_or,do_log?fplog:NULL,
1921 eprNORMAL,bCompact,mdebin,fcd,groups,&(ir->opts));
1923 if (ir->ePull != epullNO)
1925 pull_print_output(ir->pull,step,t);
1928 if (do_per_step(step,ir->nstlog))
1930 if(fflush(fplog) != 0)
1932 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of disk space?");
1938 /* Have to do this part after outputting the logfile and the edr file */
1939 state->fep_state = lamnew;
1940 for (i=0;i<efptNR;i++)
1942 state->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
1945 /* Remaining runtime */
1946 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1950 fprintf(stderr,"\n");
1952 print_time(stderr,runtime,step,ir,cr);
1955 /* Replica exchange */
1957 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1958 do_per_step(step,repl_ex_nst))
1960 bExchanged = replica_exchange(fplog,cr,repl_ex,
1964 if (bExchanged && DOMAINDECOMP(cr))
1966 dd_partition_system(fplog,step,cr,TRUE,1,
1967 state_global,top_global,ir,
1968 state,&f,mdatoms,top,fr,
1969 vsite,shellfc,constr,
1976 bStartingFromCpt = FALSE;
1978 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1979 /* With all integrators, except VV, we need to retain the pressure
1980 * at the current step for coupling at the next step.
1982 if ((state->flags & (1<<estPRES_PREV)) &&
1984 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1986 /* Store the pressure in t_state for pressure coupling
1987 * at the next MD step.
1989 copy_mat(pres,state->pres_prev);
1992 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1994 if ( (membed!=NULL) && (!bLastStep) )
1996 rescale_membed(step_rel,membed,state_global->x);
2003 /* read next frame from input trajectory */
2004 bNotLastFrame = read_next_frame(oenv,status,&rerun_fr);
2009 rerun_parallel_comm(cr,&rerun_fr,&bNotLastFrame);
2013 if (!bRerunMD || !rerun_fr.bStep)
2015 /* increase the MD step number */
2020 cycles = wallcycle_stop(wcycle,ewcSTEP);
2021 if (DOMAINDECOMP(cr) && wcycle)
2023 dd_cycles_add(cr->dd,cycles,ddCyclStep);
2026 if (bPMETuneRunning || bPMETuneTry)
2028 /* PME grid + cut-off optimization with GPUs or PME nodes */
2030 /* Count the total cycles over the last steps */
2031 cycles_pmes += cycles;
2033 /* We can only switch cut-off at NS steps */
2034 if (step % ir->nstlist == 0)
2036 /* PME grid + cut-off optimization with GPUs or PME nodes */
2039 if (DDMASTER(cr->dd))
2041 /* PME node load is too high, start tuning */
2042 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2044 dd_bcast(cr->dd,sizeof(gmx_bool),&bPMETuneRunning);
2046 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2048 bPMETuneTry = FALSE;
2051 if (bPMETuneRunning)
2053 /* init_step might not be a multiple of nstlist,
2054 * but the first cycle is always skipped anyhow.
2057 pme_load_balance(pme_loadbal,cr,
2058 (bVerbose && MASTER(cr)) ? stderr : NULL,
2060 ir,state,cycles_pmes,
2061 fr->ic,fr->nbv,&fr->pmedata,
2064 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2065 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2066 fr->rlist = fr->ic->rlist;
2067 fr->rlistlong = fr->ic->rlistlong;
2068 fr->rcoulomb = fr->ic->rcoulomb;
2069 fr->rvdw = fr->ic->rvdw;
2075 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2076 gs.set[eglsRESETCOUNTERS] != 0)
2078 /* Reset all the counters related to performance over the run */
2079 reset_all_counters(fplog,cr,step,&step_rel,ir,wcycle,nrnb,runtime,
2080 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2081 wcycle_set_reset_counters(wcycle,-1);
2082 /* Correct max_hours for the elapsed time */
2083 max_hours -= run_time/(60.0*60.0);
2084 bResetCountersHalfMaxH = FALSE;
2085 gs.set[eglsRESETCOUNTERS] = 0;
2089 /* End of main MD loop */
2093 runtime_end(runtime);
2095 if (bRerunMD && MASTER(cr))
2100 if (!(cr->duty & DUTY_PME))
2102 /* Tell the PME only node to finish */
2103 gmx_pme_send_finish(cr);
2108 if (ir->nstcalcenergy > 0 && !bRerunMD)
2110 print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
2111 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
2119 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2121 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)));
2122 fprintf(fplog,"Average number of atoms that crossed the half buffer length: %.1f\n\n",nlh.ab/nlh.nns);
2125 if (pme_loadbal != NULL)
2127 pme_loadbal_done(pme_loadbal,fplog);
2130 if (shellfc && fplog)
2132 fprintf(fplog,"Fraction of iterations that converged: %.2f %%\n",
2133 (nconverged*100.0)/step_rel);
2134 fprintf(fplog,"Average number of force evaluations per MD step: %.2f\n\n",
2138 if (repl_ex_nst > 0 && MASTER(cr))
2140 print_replica_exchange_statistics(fplog,repl_ex);
2143 runtime->nsteps_done = step_rel;