2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
58 #include "md_support.h"
59 #include "md_logging.h"
74 #include "mpelogging.h"
76 #include "domdec_network.h"
82 #include "compute_io.h"
84 #include "checkpoint.h"
85 #include "mtop_util.h"
86 #include "sighandler.h"
89 #include "pme_loadbal.h"
92 #include "types/nlistheuristics.h"
93 #include "types/iteratedconstraints.h"
94 #include "nbnxn_cuda_data_mgmt.h"
104 #include "corewrap.h"
107 static void reset_all_counters(FILE *fplog, t_commrec *cr,
108 gmx_large_int_t step,
109 gmx_large_int_t *step_rel, t_inputrec *ir,
110 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
111 gmx_runtime_t *runtime,
112 nbnxn_cuda_ptr_t cu_nbv)
114 char sbuf[STEPSTRSIZE];
116 /* Reset all the counters related to performance over the run */
117 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
118 gmx_step_str(step, sbuf));
122 nbnxn_cuda_reset_timings(cu_nbv);
125 wallcycle_stop(wcycle, ewcRUN);
126 wallcycle_reset_all(wcycle);
127 if (DOMAINDECOMP(cr))
129 reset_dd_statistics_counters(cr->dd);
132 ir->init_step += *step_rel;
133 ir->nsteps -= *step_rel;
135 wallcycle_start(wcycle, ewcRUN);
136 runtime_start(runtime);
137 print_date_and_time(fplog, cr->nodeid, "Restarted time", runtime);
140 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
141 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
143 gmx_vsite_t *vsite, gmx_constr_t constr,
144 int stepout, t_inputrec *ir,
145 gmx_mtop_t *top_global,
147 t_state *state_global,
149 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
150 gmx_edsam_t ed, t_forcerec *fr,
151 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
152 real cpt_period, real max_hours,
153 const char *deviceOptions,
155 gmx_runtime_t *runtime)
158 gmx_large_int_t step, step_rel;
160 double t, t0, lam0[efptNR];
161 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
162 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
163 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
164 bBornRadii, bStartingFromCpt;
165 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
166 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
167 bForceUpdate = FALSE, bCPT;
169 gmx_bool bMasterState;
170 int force_flags, cglo_flags;
171 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
176 t_state *bufstate = NULL;
177 matrix *scale_tot, pcoupl_mu, M, ebox;
180 gmx_repl_ex_t repl_ex = NULL;
183 t_mdebin *mdebin = NULL;
184 t_state *state = NULL;
185 rvec *f_global = NULL;
188 gmx_enerdata_t *enerd;
190 gmx_global_stat_t gstat;
191 gmx_update_t upd = NULL;
192 t_graph *graph = NULL;
194 gmx_rng_t mcrng = NULL;
196 gmx_groups_t *groups;
197 gmx_ekindata_t *ekind, *ekind_save;
198 gmx_shellfc_t shellfc;
199 int count, nconverged = 0;
202 gmx_bool bIonize = FALSE;
203 gmx_bool bTCR = FALSE, bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
205 gmx_bool bResetCountersHalfMaxH = FALSE;
206 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
207 gmx_bool bUpdateDoLR;
208 real mu_aver = 0, dvdl_constr;
209 int a0, a1, gnx = 0, ii;
210 atom_id *grpindex = NULL;
212 t_coupl_rec *tcr = NULL;
213 rvec *xcopy = NULL, *vcopy = NULL, *cbuf = NULL;
214 matrix boxcopy = {{0}}, lastbox;
216 real fom, oldfom, veta_save, pcurr, scalevir, tracevir;
223 real saved_conserved_quantity = 0;
228 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
229 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
230 gmx_iterate_t iterate;
231 gmx_large_int_t multisim_nsteps = -1; /* number of steps to do before first multisim
232 simulation stops. If equal to zero, don't
233 communicate any more between multisims.*/
234 /* PME load balancing data for GPU kernels */
235 pme_load_balancing_t pme_loadbal = NULL;
237 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
240 /* Temporary addition for FAHCORE checkpointing */
244 /* Check for special mdrun options */
245 bRerunMD = (Flags & MD_RERUN);
246 bIonize = (Flags & MD_IONIZE);
247 bFFscan = (Flags & MD_FFSCAN);
248 bAppend = (Flags & MD_APPENDFILES);
249 if (Flags & MD_RESETCOUNTERSHALFWAY)
253 /* Signal to reset the counters half the simulation steps. */
254 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
256 /* Signal to reset the counters halfway the simulation time. */
257 bResetCountersHalfMaxH = (max_hours > 0);
260 /* md-vv uses averaged full step velocities for T-control
261 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
262 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
264 if (bVV) /* to store the initial velocities while computing virial */
266 snew(cbuf, top_global->natoms);
268 /* all the iteratative cases - only if there are constraints */
269 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
270 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
271 false in this step. The correct value, true or false,
272 is set at each step, as it depends on the frequency of temperature
273 and pressure control.*/
274 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
278 /* Since we don't know if the frames read are related in any way,
279 * rebuild the neighborlist at every step.
282 ir->nstcalcenergy = 1;
286 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
288 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
289 bGStatEveryStep = (nstglobalcomm == 1);
291 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
294 "To reduce the energy communication with nstlist = -1\n"
295 "the neighbor list validity should not be checked at every step,\n"
296 "this means that exact integration is not guaranteed.\n"
297 "The neighbor list validity is checked after:\n"
298 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
299 "In most cases this will result in exact integration.\n"
300 "This reduces the energy communication by a factor of 2 to 3.\n"
301 "If you want less energy communication, set nstlist > 3.\n\n");
304 if (bRerunMD || bFFscan)
308 groups = &top_global->groups;
311 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
312 &(state_global->fep_state), lam0,
313 nrnb, top_global, &upd,
314 nfile, fnm, &outf, &mdebin,
315 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, state_global, Flags);
317 clear_mat(total_vir);
319 /* Energy terms and groups */
321 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
323 if (DOMAINDECOMP(cr))
329 snew(f, top_global->natoms);
332 /* Kinetic energy data */
334 init_ekindata(fplog, top_global, &(ir->opts), ekind);
335 /* needed for iteration of constraints */
337 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
338 /* Copy the cos acceleration to the groups struct */
339 ekind->cosacc.cos_accel = ir->cos_accel;
341 gstat = global_stat_init(ir);
344 /* Check for polarizable models and flexible constraints */
345 shellfc = init_shell_flexcon(fplog, fr->cutoff_scheme == ecutsVERLET,
346 top_global, n_flexible_constraints(constr),
347 (ir->bContinuation ||
348 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
349 NULL : state_global->x);
351 if (shellfc && ir->eI == eiNM)
353 /* Currently shells don't work with Normal Modes */
354 gmx_fatal(FARGS, "Normal Mode analysis is not supported with shells.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
357 if (vsite && ir->eI == eiNM)
359 /* Currently virtual sites don't work with Normal Modes */
360 gmx_fatal(FARGS, "Normal Mode analysis is not supported with virtual sites.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
365 #ifdef GMX_THREAD_MPI
366 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
368 set_deform_reference_box(upd,
369 deform_init_init_step_tpx,
370 deform_init_box_tpx);
371 #ifdef GMX_THREAD_MPI
372 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
377 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
378 if ((io > 2000) && MASTER(cr))
381 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
386 if (DOMAINDECOMP(cr))
388 top = dd_init_local_top(top_global);
391 dd_init_local_state(cr->dd, state_global, state);
393 if (DDMASTER(cr->dd) && ir->nstfout)
395 snew(f_global, state_global->natoms);
402 /* Initialize the particle decomposition and split the topology */
403 top = split_system(fplog, top_global, ir, cr);
405 pd_cg_range(cr, &fr->cg0, &fr->hcg);
406 pd_at_range(cr, &a0, &a1);
410 top = gmx_mtop_generate_local_top(top_global, ir);
413 a1 = top_global->natoms;
416 forcerec_set_excl_load(fr, top, cr);
418 state = partdec_init_local_state(cr, state_global);
421 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
425 set_vsite_top(vsite, top, mdatoms, cr);
428 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
430 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
435 make_local_shells(cr, mdatoms, shellfc);
438 setup_bonded_threading(fr, &top->idef);
440 if (ir->pull && PAR(cr))
442 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
446 if (DOMAINDECOMP(cr))
448 /* Distribute the charge groups over the nodes from the master node */
449 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
450 state_global, top_global, ir,
451 state, &f, mdatoms, top, fr,
452 vsite, shellfc, constr,
453 nrnb, wcycle, FALSE);
457 update_mdatoms(mdatoms, state->lambda[efptMASS]);
459 if (opt2bSet("-cpi", nfile, fnm))
461 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
465 bStateFromCP = FALSE;
470 init_expanded_ensemble(bStateFromCP,ir,&mcrng,&state->dfhist);
477 /* Update mdebin with energy history if appending to output files */
478 if (Flags & MD_APPENDFILES)
480 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
484 /* We might have read an energy history from checkpoint,
485 * free the allocated memory and reset the counts.
487 done_energyhistory(&state_global->enerhist);
488 init_energyhistory(&state_global->enerhist);
491 /* Set the initial energy history in state by updating once */
492 update_energyhistory(&state_global->enerhist, mdebin);
495 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
497 /* Set the random state if we read a checkpoint file */
498 set_stochd_state(upd, state);
501 if (state->flags & (1<<estMC_RNG))
503 set_mc_state(mcrng, state);
506 /* Initialize constraints */
509 if (!DOMAINDECOMP(cr))
511 set_constraints(constr, top, ir, mdatoms, cr);
515 /* Check whether we have to GCT stuff */
516 bTCR = ftp2bSet(efGCT, nfile, fnm);
521 fprintf(stderr, "Will do General Coupling Theory!\n");
523 gnx = top_global->mols.nr;
525 for (i = 0; (i < gnx); i++)
533 /* We need to be sure replica exchange can only occur
534 * when the energies are current */
535 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
536 "repl_ex_nst", &repl_ex_nst);
537 /* This check needs to happen before inter-simulation
538 * signals are initialized, too */
540 if (repl_ex_nst > 0 && MASTER(cr))
542 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
543 repl_ex_nst, repl_ex_nex, repl_ex_seed);
546 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
547 * With perturbed charges with soft-core we should not change the cut-off.
549 if ((Flags & MD_TUNEPME) &&
550 EEL_PME(fr->eeltype) &&
551 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
552 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
555 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
557 if (cr->duty & DUTY_PME)
559 /* Start tuning right away, as we can't measure the load */
560 bPMETuneRunning = TRUE;
564 /* Separate PME nodes, we can measure the PP/PME load balance */
569 if (!ir->bContinuation && !bRerunMD)
571 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
573 /* Set the velocities of frozen particles to zero */
574 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
576 for (m = 0; m < DIM; m++)
578 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
588 /* Constrain the initial coordinates and velocities */
589 do_constrain_first(fplog, constr, ir, mdatoms, state, f,
590 graph, cr, nrnb, fr, top, shake_vir);
594 /* Construct the virtual sites for the initial configuration */
595 construct_vsites(fplog, vsite, state->x, nrnb, ir->delta_t, NULL,
596 top->idef.iparams, top->idef.il,
597 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
603 /* set free energy calculation frequency as the minimum
604 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
605 nstfep = ir->fepvals->nstdhdl;
608 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl,nstfep);
612 nstfep = gmx_greatest_common_divisor(repl_ex_nst,nstfep);
615 /* I'm assuming we need global communication the first time! MRS */
616 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
617 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
618 | (bVV ? CGLO_PRESSURE : 0)
619 | (bVV ? CGLO_CONSTRAINT : 0)
620 | (bRerunMD ? CGLO_RERUNMD : 0)
621 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
623 bSumEkinhOld = FALSE;
624 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
625 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
626 constr, NULL, FALSE, state->box,
627 top_global, &pcurr, top_global->natoms, &bSumEkinhOld, cglo_flags);
628 if (ir->eI == eiVVAK)
630 /* a second call to get the half step temperature initialized as well */
631 /* we do the same call as above, but turn the pressure off -- internally to
632 compute_globals, this is recognized as a velocity verlet half-step
633 kinetic energy calculation. This minimized excess variables, but
634 perhaps loses some logic?*/
636 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
637 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
638 constr, NULL, FALSE, state->box,
639 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
640 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
643 /* Calculate the initial half step temperature, and save the ekinh_old */
644 if (!(Flags & MD_STARTFROMCPT))
646 for (i = 0; (i < ir->opts.ngtc); i++)
648 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
653 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
654 and there is no previous step */
657 /* if using an iterative algorithm, we need to create a working directory for the state. */
660 bufstate = init_bufstate(state);
664 snew(xcopy, state->natoms);
665 snew(vcopy, state->natoms);
666 copy_rvecn(state->x, xcopy, 0, state->natoms);
667 copy_rvecn(state->v, vcopy, 0, state->natoms);
668 copy_mat(state->box, boxcopy);
671 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
672 temperature control */
673 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
677 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
680 "RMS relative constraint deviation after constraining: %.2e\n",
681 constr_rmsd(constr, FALSE));
683 if (EI_STATE_VELOCITY(ir->eI))
685 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
689 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
690 " input trajectory '%s'\n\n",
691 *(top_global->name), opt2fn("-rerun", nfile, fnm));
694 fprintf(stderr, "Calculated time to finish depends on nsteps from "
695 "run input file,\nwhich may not correspond to the time "
696 "needed to process input trajectory.\n\n");
702 fprintf(stderr, "starting mdrun '%s'\n",
703 *(top_global->name));
706 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
710 sprintf(tbuf, "%s", "infinite");
712 if (ir->init_step > 0)
714 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
715 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
716 gmx_step_str(ir->init_step, sbuf2),
717 ir->init_step*ir->delta_t);
721 fprintf(stderr, "%s steps, %s ps.\n",
722 gmx_step_str(ir->nsteps, sbuf), tbuf);
725 fprintf(fplog, "\n");
728 print_start(fplog, cr, runtime, "mdrun");
729 runtime_start(runtime);
730 wallcycle_start(wcycle, ewcRUN);
732 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
734 chkpt_ret = fcCheckPointParallel( cr->nodeid,
738 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
743 /***********************************************************
747 ************************************************************/
749 /* if rerunMD then read coordinates and velocities from input trajectory */
752 if (getenv("GMX_FORCE_UPDATE"))
760 bNotLastFrame = read_first_frame(oenv, &status,
761 opt2fn("-rerun", nfile, fnm),
762 &rerun_fr, TRX_NEED_X | TRX_READ_V);
763 if (rerun_fr.natoms != top_global->natoms)
766 "Number of atoms in trajectory (%d) does not match the "
767 "run input file (%d)\n",
768 rerun_fr.natoms, top_global->natoms);
770 if (ir->ePBC != epbcNONE)
774 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);
776 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
778 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
785 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
788 if (ir->ePBC != epbcNONE)
790 /* Set the shift vectors.
791 * Necessary here when have a static box different from the tpr box.
793 calc_shifts(rerun_fr.box, fr->shift_vec);
797 /* loop over MD steps or if rerunMD to end of input trajectory */
799 /* Skip the first Nose-Hoover integration when we get the state from tpx */
800 bStateFromTPX = !bStateFromCP;
801 bInitStep = bFirstStep && (bStateFromTPX || bVV);
802 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
804 bSumEkinhOld = FALSE;
807 init_global_signals(&gs, cr, ir, repl_ex_nst);
809 step = ir->init_step;
812 if (ir->nstlist == -1)
814 init_nlistheuristics(&nlh, bGStatEveryStep, step);
817 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
819 /* check how many steps are left in other sims */
820 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
824 /* and stop now if we should */
825 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
826 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
827 while (!bLastStep || (bRerunMD && bNotLastFrame))
830 wallcycle_start(wcycle, ewcSTEP);
832 GMX_MPE_LOG(ev_timestep1);
838 step = rerun_fr.step;
839 step_rel = step - ir->init_step;
852 bLastStep = (step_rel == ir->nsteps);
853 t = t0 + step*ir->delta_t;
856 if (ir->efep != efepNO || ir->bSimTemp)
858 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
859 requiring different logic. */
861 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
862 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
863 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
864 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
865 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
870 update_annealing_target_temp(&(ir->opts), t);
875 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
877 for (i = 0; i < state_global->natoms; i++)
879 copy_rvec(rerun_fr.x[i], state_global->x[i]);
883 for (i = 0; i < state_global->natoms; i++)
885 copy_rvec(rerun_fr.v[i], state_global->v[i]);
890 for (i = 0; i < state_global->natoms; i++)
892 clear_rvec(state_global->v[i]);
896 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
897 " Ekin, temperature and pressure are incorrect,\n"
898 " the virial will be incorrect when constraints are present.\n"
900 bRerunWarnNoV = FALSE;
904 copy_mat(rerun_fr.box, state_global->box);
905 copy_mat(state_global->box, state->box);
907 if (vsite && (Flags & MD_RERUN_VSITE))
909 if (DOMAINDECOMP(cr))
911 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
915 /* Following is necessary because the graph may get out of sync
916 * with the coordinates if we only have every N'th coordinate set
918 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
919 shift_self(graph, state->box, state->x);
921 construct_vsites(fplog, vsite, state->x, nrnb, ir->delta_t, state->v,
922 top->idef.iparams, top->idef.il,
923 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
926 unshift_self(graph, state->box, state->x);
931 /* Stop Center of Mass motion */
932 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
934 /* Copy back starting coordinates in case we're doing a forcefield scan */
937 for (ii = 0; (ii < state->natoms); ii++)
939 copy_rvec(xcopy[ii], state->x[ii]);
940 copy_rvec(vcopy[ii], state->v[ii]);
942 copy_mat(boxcopy, state->box);
947 /* for rerun MD always do Neighbour Searching */
948 bNS = (bFirstStep || ir->nstlist != 0);
953 /* Determine whether or not to do Neighbour Searching and LR */
954 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
956 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
957 (ir->nstlist == -1 && nlh.nabnsb > 0));
959 if (bNS && ir->nstlist == -1)
961 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
965 /* check whether we should stop because another simulation has
969 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
970 (multisim_nsteps != ir->nsteps) )
977 "Stopping simulation %d because another one has finished\n",
981 gs.sig[eglsCHKPT] = 1;
986 /* < 0 means stop at next step, > 0 means stop at next NS step */
987 if ( (gs.set[eglsSTOPCOND] < 0) ||
988 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
993 /* Determine whether or not to update the Born radii if doing GB */
994 bBornRadii = bFirstStep;
995 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
1000 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
1001 do_verbose = bVerbose &&
1002 (step % stepout == 0 || bFirstStep || bLastStep);
1004 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
1008 bMasterState = TRUE;
1012 bMasterState = FALSE;
1013 /* Correct the new box if it is too skewed */
1014 if (DYNAMIC_BOX(*ir))
1016 if (correct_box(fplog, step, state->box, graph))
1018 bMasterState = TRUE;
1021 if (DOMAINDECOMP(cr) && bMasterState)
1023 dd_collect_state(cr->dd, state, state_global);
1027 if (DOMAINDECOMP(cr))
1029 /* Repartition the domain decomposition */
1030 wallcycle_start(wcycle, ewcDOMDEC);
1031 dd_partition_system(fplog, step, cr,
1032 bMasterState, nstglobalcomm,
1033 state_global, top_global, ir,
1034 state, &f, mdatoms, top, fr,
1035 vsite, shellfc, constr,
1037 do_verbose && !bPMETuneRunning);
1038 wallcycle_stop(wcycle, ewcDOMDEC);
1039 /* If using an iterative integrator, reallocate space to match the decomposition */
1043 if (MASTER(cr) && do_log && !bFFscan)
1045 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
1048 if (ir->efep != efepNO)
1050 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1053 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1056 /* We need the kinetic energy at minus the half step for determining
1057 * the full step kinetic energy and possibly for T-coupling.*/
1058 /* This may not be quite working correctly yet . . . . */
1059 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1060 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1061 constr, NULL, FALSE, state->box,
1062 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1063 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1065 clear_mat(force_vir);
1067 /* Ionize the atoms if necessary */
1070 ionize(fplog, oenv, mdatoms, top_global, t, ir, state->x, state->v,
1071 mdatoms->start, mdatoms->start+mdatoms->homenr, state->box, cr);
1074 /* Update force field in ffscan program */
1077 if (update_forcefield(fplog,
1079 mdatoms->nr, state->x, state->box))
1087 GMX_MPE_LOG(ev_timestep2);
1089 /* We write a checkpoint at this MD step when:
1090 * either at an NS step when we signalled through gs,
1091 * or at the last step (but not when we do not want confout),
1092 * but never at the first step or with rerun.
1094 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1095 (bLastStep && (Flags & MD_CONFOUT))) &&
1096 step > ir->init_step && !bRerunMD);
1099 gs.set[eglsCHKPT] = 0;
1102 /* Determine the energy and pressure:
1103 * at nstcalcenergy steps and at energy output steps (set below).
1105 if (EI_VV(ir->eI) && (!bInitStep))
1107 /* for vv, the first half of the integration actually corresponds
1108 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1109 but the virial needs to be calculated on both the current step and the 'next' step. Future
1110 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1112 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1113 bCalcVir = bCalcEner ||
1114 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1118 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1119 bCalcVir = bCalcEner ||
1120 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1123 /* Do we need global communication ? */
1124 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1125 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1126 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1128 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1130 if (do_ene || do_log)
1137 /* these CGLO_ options remain the same throughout the iteration */
1138 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1139 (bGStat ? CGLO_GSTAT : 0)
1142 force_flags = (GMX_FORCE_STATECHANGED |
1143 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1144 GMX_FORCE_ALLFORCES |
1146 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1147 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1148 (bDoFEP ? GMX_FORCE_DHDL : 0)
1153 if (do_per_step(step, ir->nstcalclr))
1155 force_flags |= GMX_FORCE_DO_LR;
1161 /* Now is the time to relax the shells */
1162 count = relax_shell_flexcon(fplog, cr, bVerbose, bFFscan ? step+1 : step,
1163 ir, bNS, force_flags,
1164 bStopCM, top, top_global,
1166 state, f, force_vir, mdatoms,
1167 nrnb, wcycle, graph, groups,
1168 shellfc, fr, bBornRadii, t, mu_tot,
1169 state->natoms, &bConverged, vsite,
1180 /* The coordinates (x) are shifted (to get whole molecules)
1182 * This is parallellized as well, and does communication too.
1183 * Check comments in sim_util.c
1185 do_force(fplog, cr, ir, step, nrnb, wcycle, top, top_global, groups,
1186 state->box, state->x, &state->hist,
1187 f, force_vir, mdatoms, enerd, fcd,
1188 state->lambda, graph,
1189 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1190 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1193 GMX_BARRIER(cr->mpi_comm_mygroup);
1197 mu_aver = calc_mu_aver(cr, state->x, mdatoms->chargeA,
1198 mu_tot, &top_global->mols, mdatoms, gnx, grpindex);
1201 if (bTCR && bFirstStep)
1203 tcr = init_coupling(fplog, nfile, fnm, cr, fr, mdatoms, &(top->idef));
1204 fprintf(fplog, "Done init_coupling\n");
1208 if (bVV && !bStartingFromCpt && !bRerunMD)
1209 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1211 if (ir->eI == eiVV && bInitStep)
1213 /* if using velocity verlet with full time step Ekin,
1214 * take the first half step only to compute the
1215 * virial for the first step. From there,
1216 * revert back to the initial coordinates
1217 * so that the input is actually the initial step.
1219 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1223 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1224 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1227 /* If we are using twin-range interactions where the long-range component
1228 * is only evaluated every nstcalclr>1 steps, we should do a special update
1229 * step to combine the long-range forces on these steps.
1230 * For nstcalclr=1 this is not done, since the forces would have been added
1231 * directly to the short-range forces already.
1233 * TODO Remove various aspects of VV+twin-range in master
1234 * branch, because VV integrators did not ever support
1235 * twin-range multiple time stepping with constraints.
1237 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1239 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1240 f, bUpdateDoLR, fr->f_twin, fcd,
1241 ekind, M, wcycle, upd, bInitStep, etrtVELOCITY1,
1242 cr, nrnb, constr, &top->idef);
1244 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1246 gmx_iterate_init(&iterate, TRUE);
1248 /* for iterations, we save these vectors, as we will be self-consistently iterating
1251 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1253 /* save the state */
1254 if (iterate.bIterationActive)
1256 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1259 bFirstIterate = TRUE;
1260 while (bFirstIterate || iterate.bIterationActive)
1262 if (iterate.bIterationActive)
1264 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1265 if (bFirstIterate && bTrotter)
1267 /* The first time through, we need a decent first estimate
1268 of veta(t+dt) to compute the constraints. Do
1269 this by computing the box volume part of the
1270 trotter integration at this time. Nothing else
1271 should be changed by this routine here. If
1272 !(first time), we start with the previous value
1275 veta_save = state->veta;
1276 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1277 vetanew = state->veta;
1278 state->veta = veta_save;
1283 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1285 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1286 state, fr->bMolPBC, graph, f,
1287 &top->idef, shake_vir, NULL,
1288 cr, nrnb, wcycle, upd, constr,
1289 bInitStep, TRUE, bCalcVir, vetanew);
1291 if (!bOK && !bFFscan)
1293 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1299 /* Need to unshift here if a do_force has been
1300 called in the previous step */
1301 unshift_self(graph, state->box, state->x);
1304 /* if VV, compute the pressure and constraints */
1305 /* For VV2, we strictly only need this if using pressure
1306 * control, but we really would like to have accurate pressures
1308 * Think about ways around this in the future?
1309 * For now, keep this choice in comments.
1311 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1312 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1314 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1315 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1317 bSumEkinhOld = TRUE;
1319 /* for vv, the first half of the integration actually corresponds to the previous step.
1320 So we need information from the last step in the first half of the integration */
1321 if (bGStat || do_per_step(step-1, nstglobalcomm))
1323 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1324 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1325 constr, NULL, FALSE, state->box,
1326 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1329 | (bTemp ? CGLO_TEMPERATURE : 0)
1330 | (bPres ? CGLO_PRESSURE : 0)
1331 | (bPres ? CGLO_CONSTRAINT : 0)
1332 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1333 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1336 /* explanation of above:
1337 a) We compute Ekin at the full time step
1338 if 1) we are using the AveVel Ekin, and it's not the
1339 initial step, or 2) if we are using AveEkin, but need the full
1340 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1341 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1342 EkinAveVel because it's needed for the pressure */
1344 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1349 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1350 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1357 /* We need the kinetic energy at minus the half step for determining
1358 * the full step kinetic energy and possibly for T-coupling.*/
1359 /* This may not be quite working correctly yet . . . . */
1360 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1361 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1362 constr, NULL, FALSE, state->box,
1363 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1364 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1369 if (iterate.bIterationActive &&
1370 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1371 state->veta, &vetanew))
1375 bFirstIterate = FALSE;
1378 if (bTrotter && !bInitStep)
1380 copy_mat(shake_vir, state->svir_prev);
1381 copy_mat(force_vir, state->fvir_prev);
1382 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1384 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1385 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE, FALSE);
1386 enerd->term[F_EKIN] = trace(ekind->ekin);
1389 /* if it's the initial step, we performed this first step just to get the constraint virial */
1390 if (bInitStep && ir->eI == eiVV)
1392 copy_rvecn(cbuf, state->v, 0, state->natoms);
1395 GMX_MPE_LOG(ev_timestep1);
1398 /* MRS -- now done iterating -- compute the conserved quantity */
1401 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1404 last_ekin = enerd->term[F_EKIN];
1406 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1408 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1410 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1413 sum_dhdl(enerd, state->lambda, ir->fepvals);
1417 /* ######## END FIRST UPDATE STEP ############## */
1418 /* ######## If doing VV, we now have v(dt) ###### */
1421 /* perform extended ensemble sampling in lambda - we don't
1422 actually move to the new state before outputting
1423 statistics, but if performing simulated tempering, we
1424 do update the velocities and the tau_t. */
1426 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, mcrng, state->v, mdatoms);
1427 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1428 copy_df_history(&state_global->dfhist,&state->dfhist);
1430 /* ################## START TRAJECTORY OUTPUT ################# */
1432 /* Now we have the energies and forces corresponding to the
1433 * coordinates at time t. We must output all of this before
1435 * for RerunMD t is read from input trajectory
1437 GMX_MPE_LOG(ev_output_start);
1440 if (do_per_step(step, ir->nstxout))
1442 mdof_flags |= MDOF_X;
1444 if (do_per_step(step, ir->nstvout))
1446 mdof_flags |= MDOF_V;
1448 if (do_per_step(step, ir->nstfout))
1450 mdof_flags |= MDOF_F;
1452 if (do_per_step(step, ir->nstxtcout))
1454 mdof_flags |= MDOF_XTC;
1458 mdof_flags |= MDOF_CPT;
1462 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1465 /* Enforce writing positions and velocities at end of run */
1466 mdof_flags |= (MDOF_X | MDOF_V);
1472 fcReportProgress( ir->nsteps, step );
1475 #if defined(__native_client__)
1476 fcCheckin(MASTER(cr));
1479 /* sync bCPT and fc record-keeping */
1480 if (bCPT && MASTER(cr))
1482 fcRequestCheckPoint();
1486 if (mdof_flags != 0)
1488 wallcycle_start(wcycle, ewcTRAJ);
1491 if (state->flags & (1<<estLD_RNG))
1493 get_stochd_state(upd, state);
1495 if (state->flags & (1<<estMC_RNG))
1497 get_mc_state(mcrng, state);
1503 state_global->ekinstate.bUpToDate = FALSE;
1507 update_ekinstate(&state_global->ekinstate, ekind);
1508 state_global->ekinstate.bUpToDate = TRUE;
1510 update_energyhistory(&state_global->enerhist, mdebin);
1513 write_traj(fplog, cr, outf, mdof_flags, top_global,
1514 step, t, state, state_global, f, f_global, &n_xtc, &x_xtc);
1521 if (bLastStep && step_rel == ir->nsteps &&
1522 (Flags & MD_CONFOUT) && MASTER(cr) &&
1523 !bRerunMD && !bFFscan)
1525 /* x and v have been collected in write_traj,
1526 * because a checkpoint file will always be written
1529 fprintf(stderr, "\nWriting final coordinates.\n");
1532 /* Make molecules whole only for confout writing */
1533 do_pbc_mtop(fplog, ir->ePBC, state->box, top_global, state_global->x);
1535 write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm),
1536 *top_global->name, top_global,
1537 state_global->x, state_global->v,
1538 ir->ePBC, state->box);
1541 wallcycle_stop(wcycle, ewcTRAJ);
1543 GMX_MPE_LOG(ev_output_finish);
1545 /* kludge -- virial is lost with restart for NPT control. Must restart */
1546 if (bStartingFromCpt && bVV)
1548 copy_mat(state->svir_prev, shake_vir);
1549 copy_mat(state->fvir_prev, force_vir);
1551 /* ################## END TRAJECTORY OUTPUT ################ */
1553 /* Determine the wallclock run time up till now */
1554 run_time = gmx_gettime() - (double)runtime->real;
1556 /* Check whether everything is still allright */
1557 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1558 #ifdef GMX_THREAD_MPI
1563 /* this is just make gs.sig compatible with the hack
1564 of sending signals around by MPI_Reduce with together with
1566 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1568 gs.sig[eglsSTOPCOND] = 1;
1570 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1572 gs.sig[eglsSTOPCOND] = -1;
1574 /* < 0 means stop at next step, > 0 means stop at next NS step */
1578 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1579 gmx_get_signal_name(),
1580 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1584 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1585 gmx_get_signal_name(),
1586 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1588 handled_stop_condition = (int)gmx_get_stop_condition();
1590 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1591 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1592 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1594 /* Signal to terminate the run */
1595 gs.sig[eglsSTOPCOND] = 1;
1598 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1600 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1603 if (bResetCountersHalfMaxH && MASTER(cr) &&
1604 run_time > max_hours*60.0*60.0*0.495)
1606 gs.sig[eglsRESETCOUNTERS] = 1;
1609 if (ir->nstlist == -1 && !bRerunMD)
1611 /* When bGStatEveryStep=FALSE, global_stat is only called
1612 * when we check the atom displacements, not at NS steps.
1613 * This means that also the bonded interaction count check is not
1614 * performed immediately after NS. Therefore a few MD steps could
1615 * be performed with missing interactions.
1616 * But wrong energies are never written to file,
1617 * since energies are only written after global_stat
1620 if (step >= nlh.step_nscheck)
1622 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1623 nlh.scale_tot, state->x);
1627 /* This is not necessarily true,
1628 * but step_nscheck is determined quite conservatively.
1634 /* In parallel we only have to check for checkpointing in steps
1635 * where we do global communication,
1636 * otherwise the other nodes don't know.
1638 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1641 run_time >= nchkpt*cpt_period*60.0)) &&
1642 gs.set[eglsCHKPT] == 0)
1644 gs.sig[eglsCHKPT] = 1;
1647 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1652 update_tcouple(fplog, step, ir, state, ekind, wcycle, upd, &MassQ, mdatoms);
1654 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1656 gmx_bool bIfRandomize;
1657 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr, DOMAINDECOMP(cr));
1658 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1659 if (constr && bIfRandomize)
1661 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1662 state, fr->bMolPBC, graph, f,
1663 &top->idef, tmp_vir, NULL,
1664 cr, nrnb, wcycle, upd, constr,
1665 bInitStep, TRUE, bCalcVir, vetanew);
1670 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1672 gmx_iterate_init(&iterate, TRUE);
1673 /* for iterations, we save these vectors, as we will be redoing the calculations */
1674 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1677 bFirstIterate = TRUE;
1678 while (bFirstIterate || iterate.bIterationActive)
1680 /* We now restore these vectors to redo the calculation with improved extended variables */
1681 if (iterate.bIterationActive)
1683 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1686 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1687 so scroll down for that logic */
1689 /* ######### START SECOND UPDATE STEP ################# */
1690 GMX_MPE_LOG(ev_update_start);
1691 /* Box is changed in update() when we do pressure coupling,
1692 * but we should still use the old box for energy corrections and when
1693 * writing it to the energy file, so it matches the trajectory files for
1694 * the same timestep above. Make a copy in a separate array.
1696 copy_mat(state->box, lastbox);
1701 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1703 wallcycle_start(wcycle, ewcUPDATE);
1704 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1707 if (iterate.bIterationActive)
1715 /* we use a new value of scalevir to converge the iterations faster */
1716 scalevir = tracevir/trace(shake_vir);
1718 msmul(shake_vir, scalevir, shake_vir);
1719 m_add(force_vir, shake_vir, total_vir);
1720 clear_mat(shake_vir);
1722 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1723 /* We can only do Berendsen coupling after we have summed
1724 * the kinetic energy or virial. Since the happens
1725 * in global_state after update, we should only do it at
1726 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1731 update_tcouple(fplog, step, ir, state, ekind, wcycle, upd, &MassQ, mdatoms);
1732 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, wcycle,
1738 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1740 /* velocity half-step update */
1741 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1742 bUpdateDoLR, fr->f_twin, fcd,
1743 ekind, M, wcycle, upd, FALSE, etrtVELOCITY2,
1744 cr, nrnb, constr, &top->idef);
1747 /* Above, initialize just copies ekinh into ekin,
1748 * it doesn't copy position (for VV),
1749 * and entire integrator for MD.
1752 if (ir->eI == eiVVAK)
1754 copy_rvecn(state->x, cbuf, 0, state->natoms);
1756 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1758 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1759 bUpdateDoLR, fr->f_twin, fcd,
1760 ekind, M, wcycle, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1761 wallcycle_stop(wcycle, ewcUPDATE);
1763 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1764 fr->bMolPBC, graph, f,
1765 &top->idef, shake_vir, force_vir,
1766 cr, nrnb, wcycle, upd, constr,
1767 bInitStep, FALSE, bCalcVir, state->veta);
1769 if (ir->eI == eiVVAK)
1771 /* erase F_EKIN and F_TEMP here? */
1772 /* just compute the kinetic energy at the half step to perform a trotter step */
1773 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1774 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1775 constr, NULL, FALSE, lastbox,
1776 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1777 cglo_flags | CGLO_TEMPERATURE
1779 wallcycle_start(wcycle, ewcUPDATE);
1780 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1781 /* now we know the scaling, we can compute the positions again again */
1782 copy_rvecn(cbuf, state->x, 0, state->natoms);
1784 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1786 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1787 bUpdateDoLR, fr->f_twin, fcd,
1788 ekind, M, wcycle, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1789 wallcycle_stop(wcycle, ewcUPDATE);
1791 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1792 /* are the small terms in the shake_vir here due
1793 * to numerical errors, or are they important
1794 * physically? I'm thinking they are just errors, but not completely sure.
1795 * For now, will call without actually constraining, constr=NULL*/
1796 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1797 state, fr->bMolPBC, graph, f,
1798 &top->idef, tmp_vir, force_vir,
1799 cr, nrnb, wcycle, upd, NULL,
1800 bInitStep, FALSE, bCalcVir,
1803 if (!bOK && !bFFscan)
1805 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1808 if (fr->bSepDVDL && fplog && do_log)
1810 fprintf(fplog, sepdvdlformat, "Constraint dV/dl", 0.0, dvdl_constr);
1814 /* this factor or 2 correction is necessary
1815 because half of the constraint force is removed
1816 in the vv step, so we have to double it. See
1817 the Redmine issue #1255. It is not yet clear
1818 if the factor of 2 is exact, or just a very
1819 good approximation, and this will be
1820 investigated. The next step is to see if this
1821 can be done adding a dhdl contribution from the
1822 rattle step, but this is somewhat more
1823 complicated with the current code. Will be
1824 investigated, hopefully for 4.6.3. However,
1825 this current solution is much better than
1826 having it completely wrong.
1828 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1832 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1837 /* Need to unshift here */
1838 unshift_self(graph, state->box, state->x);
1841 GMX_BARRIER(cr->mpi_comm_mygroup);
1842 GMX_MPE_LOG(ev_update_finish);
1846 wallcycle_start(wcycle, ewcVSITECONSTR);
1849 shift_self(graph, state->box, state->x);
1851 construct_vsites(fplog, vsite, state->x, nrnb, ir->delta_t, state->v,
1852 top->idef.iparams, top->idef.il,
1853 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1857 unshift_self(graph, state->box, state->x);
1859 wallcycle_stop(wcycle, ewcVSITECONSTR);
1862 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1863 /* With Leap-Frog we can skip compute_globals at
1864 * non-communication steps, but we need to calculate
1865 * the kinetic energy one step before communication.
1867 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1869 if (ir->nstlist == -1 && bFirstIterate)
1871 gs.sig[eglsNABNSB] = nlh.nabnsb;
1873 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1874 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1876 bFirstIterate ? &gs : NULL,
1877 (step_rel % gs.nstms == 0) &&
1878 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1880 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1882 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1883 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1884 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1885 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1886 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1887 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1890 if (ir->nstlist == -1 && bFirstIterate)
1892 nlh.nabnsb = gs.set[eglsNABNSB];
1893 gs.set[eglsNABNSB] = 0;
1896 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1897 /* ############# END CALC EKIN AND PRESSURE ################# */
1899 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1900 the virial that should probably be addressed eventually. state->veta has better properies,
1901 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1902 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1904 if (iterate.bIterationActive &&
1905 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1906 trace(shake_vir), &tracevir))
1910 bFirstIterate = FALSE;
1913 if (!bVV || bRerunMD)
1915 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1916 sum_dhdl(enerd, state->lambda, ir->fepvals);
1918 update_box(fplog, step, ir, mdatoms, state, graph, f,
1919 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, wcycle, upd, bInitStep, FALSE);
1921 /* ################# END UPDATE STEP 2 ################# */
1922 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1924 /* The coordinates (x) were unshifted in update */
1925 if (bFFscan && (shellfc == NULL || bConverged))
1927 if (print_forcefield(fplog, enerd->term, mdatoms->homenr,
1929 &(top_global->mols), mdatoms->massT, pres))
1933 fprintf(stderr, "\n");
1939 /* We will not sum ekinh_old,
1940 * so signal that we still have to do it.
1942 bSumEkinhOld = TRUE;
1947 /* Only do GCT when the relaxation of shells (minimization) has converged,
1948 * otherwise we might be coupling to bogus energies.
1949 * In parallel we must always do this, because the other sims might
1953 /* Since this is called with the new coordinates state->x, I assume
1954 * we want the new box state->box too. / EL 20040121
1956 do_coupling(fplog, oenv, nfile, fnm, tcr, t, step, enerd->term, fr,
1958 mdatoms, &(top->idef), mu_aver,
1959 top_global->mols.nr, cr,
1960 state->box, total_vir, pres,
1961 mu_tot, state->x, f, bConverged);
1965 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1967 /* use the directly determined last velocity, not actually the averaged half steps */
1968 if (bTrotter && ir->eI == eiVV)
1970 enerd->term[F_EKIN] = last_ekin;
1972 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1976 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1980 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1982 /* Check for excessively large energies */
1986 real etot_max = 1e200;
1988 real etot_max = 1e30;
1990 if (fabs(enerd->term[F_ETOT]) > etot_max)
1992 fprintf(stderr, "Energy too large (%g), giving up\n",
1993 enerd->term[F_ETOT]);
1996 /* ######### END PREPARING EDR OUTPUT ########### */
1998 /* Time for performance */
1999 if (((step % stepout) == 0) || bLastStep)
2001 runtime_upd_proc(runtime);
2007 gmx_bool do_dr, do_or;
2009 if (fplog && do_log && bDoExpanded)
2011 /* only needed if doing expanded ensemble */
2012 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
2013 &state_global->dfhist, state->fep_state, ir->nstlog, step);
2015 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
2019 upd_mdebin(mdebin, bDoDHDL, TRUE,
2020 t, mdatoms->tmass, enerd, state,
2021 ir->fepvals, ir->expandedvals, lastbox,
2022 shake_vir, force_vir, total_vir, pres,
2023 ekind, mu_tot, constr);
2027 upd_mdebin_step(mdebin);
2030 do_dr = do_per_step(step, ir->nstdisreout);
2031 do_or = do_per_step(step, ir->nstorireout);
2033 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
2035 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
2037 if (ir->ePull != epullNO)
2039 pull_print_output(ir->pull, step, t);
2042 if (do_per_step(step, ir->nstlog))
2044 if (fflush(fplog) != 0)
2046 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
2052 /* Have to do this part _after_ outputting the logfile and the edr file */
2053 /* Gets written into the state at the beginning of next loop*/
2054 state->fep_state = lamnew;
2057 /* Remaining runtime */
2058 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
2062 fprintf(stderr, "\n");
2064 print_time(stderr, runtime, step, ir, cr);
2067 /* Replica exchange */
2069 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
2070 do_per_step(step, repl_ex_nst))
2072 bExchanged = replica_exchange(fplog, cr, repl_ex,
2073 state_global, enerd,
2076 if (bExchanged && DOMAINDECOMP(cr))
2078 dd_partition_system(fplog, step, cr, TRUE, 1,
2079 state_global, top_global, ir,
2080 state, &f, mdatoms, top, fr,
2081 vsite, shellfc, constr,
2082 nrnb, wcycle, FALSE);
2088 bStartingFromCpt = FALSE;
2090 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
2091 /* With all integrators, except VV, we need to retain the pressure
2092 * at the current step for coupling at the next step.
2094 if ((state->flags & (1<<estPRES_PREV)) &&
2096 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
2098 /* Store the pressure in t_state for pressure coupling
2099 * at the next MD step.
2101 copy_mat(pres, state->pres_prev);
2104 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2106 if ( (membed != NULL) && (!bLastStep) )
2108 rescale_membed(step_rel, membed, state_global->x);
2115 /* read next frame from input trajectory */
2116 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
2121 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
2125 if (!bRerunMD || !rerun_fr.bStep)
2127 /* increase the MD step number */
2132 cycles = wallcycle_stop(wcycle, ewcSTEP);
2133 if (DOMAINDECOMP(cr) && wcycle)
2135 dd_cycles_add(cr->dd, cycles, ddCyclStep);
2138 if (bPMETuneRunning || bPMETuneTry)
2140 /* PME grid + cut-off optimization with GPUs or PME nodes */
2142 /* Count the total cycles over the last steps */
2143 cycles_pmes += cycles;
2145 /* We can only switch cut-off at NS steps */
2146 if (step % ir->nstlist == 0)
2148 /* PME grid + cut-off optimization with GPUs or PME nodes */
2151 if (DDMASTER(cr->dd))
2153 /* PME node load is too high, start tuning */
2154 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2156 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
2158 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2160 bPMETuneTry = FALSE;
2163 if (bPMETuneRunning)
2165 /* init_step might not be a multiple of nstlist,
2166 * but the first cycle is always skipped anyhow.
2169 pme_load_balance(pme_loadbal, cr,
2170 (bVerbose && MASTER(cr)) ? stderr : NULL,
2172 ir, state, cycles_pmes,
2173 fr->ic, fr->nbv, &fr->pmedata,
2176 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2177 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2178 fr->rlist = fr->ic->rlist;
2179 fr->rlistlong = fr->ic->rlistlong;
2180 fr->rcoulomb = fr->ic->rcoulomb;
2181 fr->rvdw = fr->ic->rvdw;
2187 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2188 gs.set[eglsRESETCOUNTERS] != 0)
2190 /* Reset all the counters related to performance over the run */
2191 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, runtime,
2192 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2193 wcycle_set_reset_counters(wcycle, -1);
2194 if (!(cr->duty & DUTY_PME))
2196 /* Tell our PME node to reset its counters */
2197 gmx_pme_send_resetcounters(cr, step);
2199 /* Correct max_hours for the elapsed time */
2200 max_hours -= run_time/(60.0*60.0);
2201 bResetCountersHalfMaxH = FALSE;
2202 gs.set[eglsRESETCOUNTERS] = 0;
2206 /* End of main MD loop */
2210 runtime_end(runtime);
2212 if (bRerunMD && MASTER(cr))
2217 if (!(cr->duty & DUTY_PME))
2219 /* Tell the PME only node to finish */
2220 gmx_pme_send_finish(cr);
2225 if (ir->nstcalcenergy > 0 && !bRerunMD)
2227 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
2228 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
2236 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2238 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)));
2239 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
2242 if (pme_loadbal != NULL)
2244 pme_loadbal_done(pme_loadbal, cr, fplog,
2245 fr->nbv != NULL && fr->nbv->bUseGPU);
2248 if (shellfc && fplog)
2250 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
2251 (nconverged*100.0)/step_rel);
2252 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
2256 if (repl_ex_nst > 0 && MASTER(cr))
2258 print_replica_exchange_statistics(fplog, repl_ex);
2261 runtime->nsteps_done = step_rel;