1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
56 #include "md_support.h"
57 #include "md_logging.h"
73 #include "domdec_network.h"
79 #include "compute_io.h"
81 #include "checkpoint.h"
82 #include "mtop_util.h"
83 #include "sighandler.h"
86 #include "pme_loadbal.h"
89 #include "types/nlistheuristics.h"
90 #include "types/iteratedconstraints.h"
91 #include "nbnxn_cuda_data_mgmt.h"
101 #include "corewrap.h"
104 static void reset_all_counters(FILE *fplog, t_commrec *cr,
105 gmx_large_int_t step,
106 gmx_large_int_t *step_rel, t_inputrec *ir,
107 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
108 gmx_runtime_t *runtime,
109 nbnxn_cuda_ptr_t cu_nbv)
111 char sbuf[STEPSTRSIZE];
113 /* Reset all the counters related to performance over the run */
114 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
115 gmx_step_str(step, sbuf));
119 nbnxn_cuda_reset_timings(cu_nbv);
122 wallcycle_stop(wcycle, ewcRUN);
123 wallcycle_reset_all(wcycle);
124 if (DOMAINDECOMP(cr))
126 reset_dd_statistics_counters(cr->dd);
129 ir->init_step += *step_rel;
130 ir->nsteps -= *step_rel;
132 wallcycle_start(wcycle, ewcRUN);
133 runtime_start(runtime);
134 print_date_and_time(fplog, cr->nodeid, "Restarted time", runtime);
137 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
138 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
140 gmx_vsite_t *vsite, gmx_constr_t constr,
141 int stepout, t_inputrec *ir,
142 gmx_mtop_t *top_global,
144 t_state *state_global,
146 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
147 gmx_edsam_t ed, t_forcerec *fr,
148 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
149 real cpt_period, real max_hours,
150 const char *deviceOptions,
152 gmx_runtime_t *runtime)
155 gmx_large_int_t step, step_rel;
157 double t, t0, lam0[efptNR];
158 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
159 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
160 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
161 bBornRadii, bStartingFromCpt;
162 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
163 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
164 bForceUpdate = FALSE, bCPT;
166 gmx_bool bMasterState;
167 int force_flags, cglo_flags;
168 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
173 t_state *bufstate = NULL;
174 matrix *scale_tot, pcoupl_mu, M, ebox;
177 gmx_repl_ex_t repl_ex = NULL;
180 t_mdebin *mdebin = NULL;
181 df_history_t df_history;
182 t_state *state = NULL;
183 rvec *f_global = NULL;
186 gmx_enerdata_t *enerd;
188 gmx_global_stat_t gstat;
189 gmx_update_t upd = NULL;
190 t_graph *graph = NULL;
192 gmx_rng_t mcrng = NULL;
194 gmx_groups_t *groups;
195 gmx_ekindata_t *ekind, *ekind_save;
196 gmx_shellfc_t shellfc;
197 int count, nconverged = 0;
200 gmx_bool bIonize = FALSE;
201 gmx_bool bTCR = FALSE, bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
203 gmx_bool bResetCountersHalfMaxH = FALSE;
204 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
205 gmx_bool bUpdateDoLR;
206 real mu_aver = 0, dvdl;
207 int a0, a1, gnx = 0, ii;
208 atom_id *grpindex = NULL;
210 t_coupl_rec *tcr = NULL;
211 rvec *xcopy = NULL, *vcopy = NULL, *cbuf = NULL;
212 matrix boxcopy = {{0}}, lastbox;
214 real fom, oldfom, veta_save, pcurr, scalevir, tracevir;
221 real saved_conserved_quantity = 0;
226 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
227 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
228 gmx_iterate_t iterate;
229 gmx_large_int_t multisim_nsteps = -1; /* number of steps to do before first multisim
230 simulation stops. If equal to zero, don't
231 communicate any more between multisims.*/
232 /* PME load balancing data for GPU kernels */
233 pme_load_balancing_t pme_loadbal = NULL;
235 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
238 /* Temporary addition for FAHCORE checkpointing */
242 /* Check for special mdrun options */
243 bRerunMD = (Flags & MD_RERUN);
244 bIonize = (Flags & MD_IONIZE);
245 bFFscan = (Flags & MD_FFSCAN);
246 bAppend = (Flags & MD_APPENDFILES);
247 if (Flags & MD_RESETCOUNTERSHALFWAY)
251 /* Signal to reset the counters half the simulation steps. */
252 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
254 /* Signal to reset the counters halfway the simulation time. */
255 bResetCountersHalfMaxH = (max_hours > 0);
258 /* md-vv uses averaged full step velocities for T-control
259 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
260 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
262 if (bVV) /* to store the initial velocities while computing virial */
264 snew(cbuf, top_global->natoms);
266 /* all the iteratative cases - only if there are constraints */
267 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
268 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
269 false in this step. The correct value, true or false,
270 is set at each step, as it depends on the frequency of temperature
271 and pressure control.*/
272 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
276 /* Since we don't know if the frames read are related in any way,
277 * rebuild the neighborlist at every step.
280 ir->nstcalcenergy = 1;
284 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
286 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
287 bGStatEveryStep = (nstglobalcomm == 1);
289 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
292 "To reduce the energy communication with nstlist = -1\n"
293 "the neighbor list validity should not be checked at every step,\n"
294 "this means that exact integration is not guaranteed.\n"
295 "The neighbor list validity is checked after:\n"
296 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
297 "In most cases this will result in exact integration.\n"
298 "This reduces the energy communication by a factor of 2 to 3.\n"
299 "If you want less energy communication, set nstlist > 3.\n\n");
302 if (bRerunMD || bFFscan)
306 groups = &top_global->groups;
309 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
310 &(state_global->fep_state), lam0,
311 nrnb, top_global, &upd,
312 nfile, fnm, &outf, &mdebin,
313 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, state_global, Flags);
315 clear_mat(total_vir);
317 /* Energy terms and groups */
319 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
321 if (DOMAINDECOMP(cr))
327 snew(f, top_global->natoms);
330 /* lambda Monte carlo random number generator */
333 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
335 /* copy the state into df_history */
336 copy_df_history(&df_history, &state_global->dfhist);
338 /* Kinetic energy data */
340 init_ekindata(fplog, top_global, &(ir->opts), ekind);
341 /* needed for iteration of constraints */
343 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
344 /* Copy the cos acceleration to the groups struct */
345 ekind->cosacc.cos_accel = ir->cos_accel;
347 gstat = global_stat_init(ir);
350 /* Check for polarizable models and flexible constraints */
351 shellfc = init_shell_flexcon(fplog,
352 top_global, n_flexible_constraints(constr),
353 (ir->bContinuation ||
354 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
355 NULL : state_global->x);
359 #ifdef GMX_THREAD_MPI
360 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
362 set_deform_reference_box(upd,
363 deform_init_init_step_tpx,
364 deform_init_box_tpx);
365 #ifdef GMX_THREAD_MPI
366 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
371 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
372 if ((io > 2000) && MASTER(cr))
375 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
380 if (DOMAINDECOMP(cr))
382 top = dd_init_local_top(top_global);
385 dd_init_local_state(cr->dd, state_global, state);
387 if (DDMASTER(cr->dd) && ir->nstfout)
389 snew(f_global, state_global->natoms);
396 /* Initialize the particle decomposition and split the topology */
397 top = split_system(fplog, top_global, ir, cr);
399 pd_cg_range(cr, &fr->cg0, &fr->hcg);
400 pd_at_range(cr, &a0, &a1);
404 top = gmx_mtop_generate_local_top(top_global, ir);
407 a1 = top_global->natoms;
410 forcerec_set_excl_load(fr, top, cr);
412 state = partdec_init_local_state(cr, state_global);
415 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
419 set_vsite_top(vsite, top, mdatoms, cr);
422 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
424 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
429 make_local_shells(cr, mdatoms, shellfc);
432 init_bonded_thread_force_reduction(fr, &top->idef);
434 if (ir->pull && PAR(cr))
436 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
440 if (DOMAINDECOMP(cr))
442 /* Distribute the charge groups over the nodes from the master node */
443 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
444 state_global, top_global, ir,
445 state, &f, mdatoms, top, fr,
446 vsite, shellfc, constr,
447 nrnb, wcycle, FALSE);
451 update_mdatoms(mdatoms, state->lambda[efptMASS]);
453 if (opt2bSet("-cpi", nfile, fnm))
455 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
459 bStateFromCP = FALSE;
466 /* Update mdebin with energy history if appending to output files */
467 if (Flags & MD_APPENDFILES)
469 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
473 /* We might have read an energy history from checkpoint,
474 * free the allocated memory and reset the counts.
476 done_energyhistory(&state_global->enerhist);
477 init_energyhistory(&state_global->enerhist);
480 /* Set the initial energy history in state by updating once */
481 update_energyhistory(&state_global->enerhist, mdebin);
484 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
486 /* Set the random state if we read a checkpoint file */
487 set_stochd_state(upd, state);
490 if (state->flags & (1<<estMC_RNG))
492 set_mc_state(mcrng, state);
495 /* Initialize constraints */
498 if (!DOMAINDECOMP(cr))
500 set_constraints(constr, top, ir, mdatoms, cr);
504 /* Check whether we have to GCT stuff */
505 bTCR = ftp2bSet(efGCT, nfile, fnm);
510 fprintf(stderr, "Will do General Coupling Theory!\n");
512 gnx = top_global->mols.nr;
514 for (i = 0; (i < gnx); i++)
522 /* We need to be sure replica exchange can only occur
523 * when the energies are current */
524 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
525 "repl_ex_nst", &repl_ex_nst);
526 /* This check needs to happen before inter-simulation
527 * signals are initialized, too */
529 if (repl_ex_nst > 0 && MASTER(cr))
531 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
532 repl_ex_nst, repl_ex_nex, repl_ex_seed);
535 /* PME tuning is only supported with GPUs or PME nodes and not with rerun */
536 if ((Flags & MD_TUNEPME) &&
537 EEL_PME(fr->eeltype) &&
538 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
541 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
543 if (cr->duty & DUTY_PME)
545 /* Start tuning right away, as we can't measure the load */
546 bPMETuneRunning = TRUE;
550 /* Separate PME nodes, we can measure the PP/PME load balance */
555 if (!ir->bContinuation && !bRerunMD)
557 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
559 /* Set the velocities of frozen particles to zero */
560 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
562 for (m = 0; m < DIM; m++)
564 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
574 /* Constrain the initial coordinates and velocities */
575 do_constrain_first(fplog, constr, ir, mdatoms, state, f,
576 graph, cr, nrnb, fr, top, shake_vir);
580 /* Construct the virtual sites for the initial configuration */
581 construct_vsites(fplog, vsite, state->x, nrnb, ir->delta_t, NULL,
582 top->idef.iparams, top->idef.il,
583 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
589 /* set free energy calculation frequency as the minimum of nstdhdl, nstexpanded, and nstrepl_ex_nst*/
590 nstfep = ir->fepvals->nstdhdl;
591 if (ir->bExpanded && (nstfep > ir->expandedvals->nstexpanded))
593 nstfep = ir->expandedvals->nstexpanded;
595 if (repl_ex_nst > 0 && nstfep > repl_ex_nst)
597 nstfep = repl_ex_nst;
600 /* I'm assuming we need global communication the first time! MRS */
601 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
602 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
603 | (bVV ? CGLO_PRESSURE : 0)
604 | (bVV ? CGLO_CONSTRAINT : 0)
605 | (bRerunMD ? CGLO_RERUNMD : 0)
606 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
608 bSumEkinhOld = FALSE;
609 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
610 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
611 constr, NULL, FALSE, state->box,
612 top_global, &pcurr, top_global->natoms, &bSumEkinhOld, cglo_flags);
613 if (ir->eI == eiVVAK)
615 /* a second call to get the half step temperature initialized as well */
616 /* we do the same call as above, but turn the pressure off -- internally to
617 compute_globals, this is recognized as a velocity verlet half-step
618 kinetic energy calculation. This minimized excess variables, but
619 perhaps loses some logic?*/
621 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
622 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
623 constr, NULL, FALSE, state->box,
624 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
625 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
628 /* Calculate the initial half step temperature, and save the ekinh_old */
629 if (!(Flags & MD_STARTFROMCPT))
631 for (i = 0; (i < ir->opts.ngtc); i++)
633 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
638 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
639 and there is no previous step */
642 /* if using an iterative algorithm, we need to create a working directory for the state. */
645 bufstate = init_bufstate(state);
649 snew(xcopy, state->natoms);
650 snew(vcopy, state->natoms);
651 copy_rvecn(state->x, xcopy, 0, state->natoms);
652 copy_rvecn(state->v, vcopy, 0, state->natoms);
653 copy_mat(state->box, boxcopy);
656 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
657 temperature control */
658 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
662 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
665 "RMS relative constraint deviation after constraining: %.2e\n",
666 constr_rmsd(constr, FALSE));
668 if (EI_STATE_VELOCITY(ir->eI))
670 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
674 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
675 " input trajectory '%s'\n\n",
676 *(top_global->name), opt2fn("-rerun", nfile, fnm));
679 fprintf(stderr, "Calculated time to finish depends on nsteps from "
680 "run input file,\nwhich may not correspond to the time "
681 "needed to process input trajectory.\n\n");
687 fprintf(stderr, "starting mdrun '%s'\n",
688 *(top_global->name));
691 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
695 sprintf(tbuf, "%s", "infinite");
697 if (ir->init_step > 0)
699 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
700 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
701 gmx_step_str(ir->init_step, sbuf2),
702 ir->init_step*ir->delta_t);
706 fprintf(stderr, "%s steps, %s ps.\n",
707 gmx_step_str(ir->nsteps, sbuf), tbuf);
710 fprintf(fplog, "\n");
713 /* Set and write start time */
714 runtime_start(runtime);
715 print_date_and_time(fplog, cr->nodeid, "Started mdrun", runtime);
716 wallcycle_start(wcycle, ewcRUN);
719 fprintf(fplog, "\n");
722 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
724 chkpt_ret = fcCheckPointParallel( cr->nodeid,
728 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
733 /***********************************************************
737 ************************************************************/
739 /* if rerunMD then read coordinates and velocities from input trajectory */
742 if (getenv("GMX_FORCE_UPDATE"))
750 bNotLastFrame = read_first_frame(oenv, &status,
751 opt2fn("-rerun", nfile, fnm),
752 &rerun_fr, TRX_NEED_X | TRX_READ_V);
753 if (rerun_fr.natoms != top_global->natoms)
756 "Number of atoms in trajectory (%d) does not match the "
757 "run input file (%d)\n",
758 rerun_fr.natoms, top_global->natoms);
760 if (ir->ePBC != epbcNONE)
764 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);
766 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
768 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
775 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
778 if (ir->ePBC != epbcNONE)
780 /* Set the shift vectors.
781 * Necessary here when have a static box different from the tpr box.
783 calc_shifts(rerun_fr.box, fr->shift_vec);
787 /* loop over MD steps or if rerunMD to end of input trajectory */
789 /* Skip the first Nose-Hoover integration when we get the state from tpx */
790 bStateFromTPX = !bStateFromCP;
791 bInitStep = bFirstStep && (bStateFromTPX || bVV);
792 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
794 bSumEkinhOld = FALSE;
797 init_global_signals(&gs, cr, ir, repl_ex_nst);
799 step = ir->init_step;
802 if (ir->nstlist == -1)
804 init_nlistheuristics(&nlh, bGStatEveryStep, step);
807 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
809 /* check how many steps are left in other sims */
810 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
814 /* and stop now if we should */
815 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
816 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
817 while (!bLastStep || (bRerunMD && bNotLastFrame))
820 wallcycle_start(wcycle, ewcSTEP);
826 step = rerun_fr.step;
827 step_rel = step - ir->init_step;
840 bLastStep = (step_rel == ir->nsteps);
841 t = t0 + step*ir->delta_t;
844 if (ir->efep != efepNO || ir->bSimTemp)
846 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
847 requiring different logic. */
849 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
850 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
851 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
852 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
857 update_annealing_target_temp(&(ir->opts), t);
862 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
864 for (i = 0; i < state_global->natoms; i++)
866 copy_rvec(rerun_fr.x[i], state_global->x[i]);
870 for (i = 0; i < state_global->natoms; i++)
872 copy_rvec(rerun_fr.v[i], state_global->v[i]);
877 for (i = 0; i < state_global->natoms; i++)
879 clear_rvec(state_global->v[i]);
883 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
884 " Ekin, temperature and pressure are incorrect,\n"
885 " the virial will be incorrect when constraints are present.\n"
887 bRerunWarnNoV = FALSE;
891 copy_mat(rerun_fr.box, state_global->box);
892 copy_mat(state_global->box, state->box);
894 if (vsite && (Flags & MD_RERUN_VSITE))
896 if (DOMAINDECOMP(cr))
898 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
902 /* Following is necessary because the graph may get out of sync
903 * with the coordinates if we only have every N'th coordinate set
905 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
906 shift_self(graph, state->box, state->x);
908 construct_vsites(fplog, vsite, state->x, nrnb, ir->delta_t, state->v,
909 top->idef.iparams, top->idef.il,
910 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
913 unshift_self(graph, state->box, state->x);
918 /* Stop Center of Mass motion */
919 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
921 /* Copy back starting coordinates in case we're doing a forcefield scan */
924 for (ii = 0; (ii < state->natoms); ii++)
926 copy_rvec(xcopy[ii], state->x[ii]);
927 copy_rvec(vcopy[ii], state->v[ii]);
929 copy_mat(boxcopy, state->box);
934 /* for rerun MD always do Neighbour Searching */
935 bNS = (bFirstStep || ir->nstlist != 0);
940 /* Determine whether or not to do Neighbour Searching and LR */
941 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
943 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
944 (ir->nstlist == -1 && nlh.nabnsb > 0));
946 if (bNS && ir->nstlist == -1)
948 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
952 /* check whether we should stop because another simulation has
956 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
957 (multisim_nsteps != ir->nsteps) )
964 "Stopping simulation %d because another one has finished\n",
968 gs.sig[eglsCHKPT] = 1;
973 /* < 0 means stop at next step, > 0 means stop at next NS step */
974 if ( (gs.set[eglsSTOPCOND] < 0 ) ||
975 ( (gs.set[eglsSTOPCOND] > 0 ) && ( bNS || ir->nstlist == 0)) )
980 /* Determine whether or not to update the Born radii if doing GB */
981 bBornRadii = bFirstStep;
982 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
987 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
988 do_verbose = bVerbose &&
989 (step % stepout == 0 || bFirstStep || bLastStep);
991 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
999 bMasterState = FALSE;
1000 /* Correct the new box if it is too skewed */
1001 if (DYNAMIC_BOX(*ir))
1003 if (correct_box(fplog, step, state->box, graph))
1005 bMasterState = TRUE;
1008 if (DOMAINDECOMP(cr) && bMasterState)
1010 dd_collect_state(cr->dd, state, state_global);
1014 if (DOMAINDECOMP(cr))
1016 /* Repartition the domain decomposition */
1017 wallcycle_start(wcycle, ewcDOMDEC);
1018 dd_partition_system(fplog, step, cr,
1019 bMasterState, nstglobalcomm,
1020 state_global, top_global, ir,
1021 state, &f, mdatoms, top, fr,
1022 vsite, shellfc, constr,
1024 do_verbose && !bPMETuneRunning);
1025 wallcycle_stop(wcycle, ewcDOMDEC);
1026 /* If using an iterative integrator, reallocate space to match the decomposition */
1030 if (MASTER(cr) && do_log && !bFFscan)
1032 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
1035 if (ir->efep != efepNO)
1037 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1040 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1043 /* We need the kinetic energy at minus the half step for determining
1044 * the full step kinetic energy and possibly for T-coupling.*/
1045 /* This may not be quite working correctly yet . . . . */
1046 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1047 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1048 constr, NULL, FALSE, state->box,
1049 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1050 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1052 clear_mat(force_vir);
1054 /* Ionize the atoms if necessary */
1057 ionize(fplog, oenv, mdatoms, top_global, t, ir, state->x, state->v,
1058 mdatoms->start, mdatoms->start+mdatoms->homenr, state->box, cr);
1061 /* Update force field in ffscan program */
1064 if (update_forcefield(fplog,
1066 mdatoms->nr, state->x, state->box))
1074 /* We write a checkpoint at this MD step when:
1075 * either at an NS step when we signalled through gs,
1076 * or at the last step (but not when we do not want confout),
1077 * but never at the first step or with rerun.
1079 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1080 (bLastStep && (Flags & MD_CONFOUT))) &&
1081 step > ir->init_step && !bRerunMD);
1084 gs.set[eglsCHKPT] = 0;
1087 /* Determine the energy and pressure:
1088 * at nstcalcenergy steps and at energy output steps (set below).
1090 if (EI_VV(ir->eI) && (!bInitStep))
1092 /* for vv, the first half of the integration actually corresponds
1093 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1094 but the virial needs to be calculated on both the current step and the 'next' step. Future
1095 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1097 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1098 bCalcVir = bCalcEner ||
1099 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1103 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1104 bCalcVir = bCalcEner ||
1105 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1108 /* Do we need global communication ? */
1109 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1110 do_per_step(step, nstglobalcomm) ||
1111 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1113 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1115 if (do_ene || do_log)
1122 /* these CGLO_ options remain the same throughout the iteration */
1123 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1124 (bGStat ? CGLO_GSTAT : 0)
1127 force_flags = (GMX_FORCE_STATECHANGED |
1128 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1129 GMX_FORCE_ALLFORCES |
1131 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1132 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1133 (bDoFEP ? GMX_FORCE_DHDL : 0)
1138 if (do_per_step(step, ir->nstcalclr))
1140 force_flags |= GMX_FORCE_DO_LR;
1146 /* Now is the time to relax the shells */
1147 count = relax_shell_flexcon(fplog, cr, bVerbose, bFFscan ? step+1 : step,
1148 ir, bNS, force_flags,
1149 bStopCM, top, top_global,
1151 state, f, force_vir, mdatoms,
1152 nrnb, wcycle, graph, groups,
1153 shellfc, fr, bBornRadii, t, mu_tot,
1154 state->natoms, &bConverged, vsite,
1165 /* The coordinates (x) are shifted (to get whole molecules)
1167 * This is parallellized as well, and does communication too.
1168 * Check comments in sim_util.c
1170 do_force(fplog, cr, ir, step, nrnb, wcycle, top, top_global, groups,
1171 state->box, state->x, &state->hist,
1172 f, force_vir, mdatoms, enerd, fcd,
1173 state->lambda, graph,
1174 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1175 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1180 mu_aver = calc_mu_aver(cr, state->x, mdatoms->chargeA,
1181 mu_tot, &top_global->mols, mdatoms, gnx, grpindex);
1184 if (bTCR && bFirstStep)
1186 tcr = init_coupling(fplog, nfile, fnm, cr, fr, mdatoms, &(top->idef));
1187 fprintf(fplog, "Done init_coupling\n");
1191 if (bVV && !bStartingFromCpt && !bRerunMD)
1192 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1194 if (ir->eI == eiVV && bInitStep)
1196 /* if using velocity verlet with full time step Ekin,
1197 * take the first half step only to compute the
1198 * virial for the first step. From there,
1199 * revert back to the initial coordinates
1200 * so that the input is actually the initial step.
1202 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1206 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1207 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1210 /* If we are using twin-range interactions where the long-range component
1211 * is only evaluated every nstcalclr>1 steps, we should do a special update
1212 * step to combine the long-range forces on these steps.
1213 * For nstcalclr=1 this is not done, since the forces would have been added
1214 * directly to the short-range forces already.
1216 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1218 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1219 f, bUpdateDoLR, fr->f_twin, fcd,
1220 ekind, M, wcycle, upd, bInitStep, etrtVELOCITY1,
1221 cr, nrnb, constr, &top->idef);
1223 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1225 gmx_iterate_init(&iterate, TRUE);
1227 /* for iterations, we save these vectors, as we will be self-consistently iterating
1230 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1232 /* save the state */
1233 if (iterate.bIterationActive)
1235 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1238 bFirstIterate = TRUE;
1239 while (bFirstIterate || iterate.bIterationActive)
1241 if (iterate.bIterationActive)
1243 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1244 if (bFirstIterate && bTrotter)
1246 /* The first time through, we need a decent first estimate
1247 of veta(t+dt) to compute the constraints. Do
1248 this by computing the box volume part of the
1249 trotter integration at this time. Nothing else
1250 should be changed by this routine here. If
1251 !(first time), we start with the previous value
1254 veta_save = state->veta;
1255 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1256 vetanew = state->veta;
1257 state->veta = veta_save;
1262 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1266 update_constraints(fplog, step, &dvdl, ir, ekind, mdatoms,
1267 state, fr->bMolPBC, graph, f,
1268 &top->idef, shake_vir, NULL,
1269 cr, nrnb, wcycle, upd, constr,
1270 bInitStep, TRUE, bCalcVir, vetanew);
1272 if (!bOK && !bFFscan)
1274 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1280 /* Need to unshift here if a do_force has been
1281 called in the previous step */
1282 unshift_self(graph, state->box, state->x);
1285 /* if VV, compute the pressure and constraints */
1286 /* For VV2, we strictly only need this if using pressure
1287 * control, but we really would like to have accurate pressures
1289 * Think about ways around this in the future?
1290 * For now, keep this choice in comments.
1292 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1293 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1295 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1296 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1298 bSumEkinhOld = TRUE;
1300 /* for vv, the first half of the integration actually corresponds to the previous step.
1301 So we need information from the last step in the first half of the integration */
1302 if (bGStat || do_per_step(step-1, nstglobalcomm))
1304 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1305 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1306 constr, NULL, FALSE, state->box,
1307 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1310 | (bTemp ? CGLO_TEMPERATURE : 0)
1311 | (bPres ? CGLO_PRESSURE : 0)
1312 | (bPres ? CGLO_CONSTRAINT : 0)
1313 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1314 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1317 /* explanation of above:
1318 a) We compute Ekin at the full time step
1319 if 1) we are using the AveVel Ekin, and it's not the
1320 initial step, or 2) if we are using AveEkin, but need the full
1321 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1322 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1323 EkinAveVel because it's needed for the pressure */
1325 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1330 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1331 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1338 /* We need the kinetic energy at minus the half step for determining
1339 * the full step kinetic energy and possibly for T-coupling.*/
1340 /* This may not be quite working correctly yet . . . . */
1341 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1342 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1343 constr, NULL, FALSE, state->box,
1344 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1345 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1350 if (iterate.bIterationActive &&
1351 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1352 state->veta, &vetanew))
1356 bFirstIterate = FALSE;
1359 if (bTrotter && !bInitStep)
1361 enerd->term[F_DVDL_BONDED] += dvdl; /* only add after iterations */
1362 copy_mat(shake_vir, state->svir_prev);
1363 copy_mat(force_vir, state->fvir_prev);
1364 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1366 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1367 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE, FALSE);
1368 enerd->term[F_EKIN] = trace(ekind->ekin);
1371 /* if it's the initial step, we performed this first step just to get the constraint virial */
1372 if (bInitStep && ir->eI == eiVV)
1374 copy_rvecn(cbuf, state->v, 0, state->natoms);
1377 if (fr->bSepDVDL && fplog && do_log)
1379 fprintf(fplog, sepdvdlformat, "Constraint", 0.0, dvdl);
1381 enerd->term[F_DVDL_BONDED] += dvdl;
1384 /* MRS -- now done iterating -- compute the conserved quantity */
1387 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1390 last_ekin = enerd->term[F_EKIN];
1392 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1394 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1396 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1399 sum_dhdl(enerd, state->lambda, ir->fepvals);
1403 /* ######## END FIRST UPDATE STEP ############## */
1404 /* ######## If doing VV, we now have v(dt) ###### */
1407 /* perform extended ensemble sampling in lambda - we don't
1408 actually move to the new state before outputting
1409 statistics, but if performing simulated tempering, we
1410 do update the velocities and the tau_t. */
1412 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, &df_history, step, mcrng, state->v, mdatoms);
1414 /* ################## START TRAJECTORY OUTPUT ################# */
1416 /* Now we have the energies and forces corresponding to the
1417 * coordinates at time t. We must output all of this before
1419 * for RerunMD t is read from input trajectory
1422 if (do_per_step(step, ir->nstxout))
1424 mdof_flags |= MDOF_X;
1426 if (do_per_step(step, ir->nstvout))
1428 mdof_flags |= MDOF_V;
1430 if (do_per_step(step, ir->nstfout))
1432 mdof_flags |= MDOF_F;
1434 if (do_per_step(step, ir->nstxtcout))
1436 mdof_flags |= MDOF_XTC;
1440 mdof_flags |= MDOF_CPT;
1444 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1447 /* Enforce writing positions and velocities at end of run */
1448 mdof_flags |= (MDOF_X | MDOF_V);
1454 fcReportProgress( ir->nsteps, step );
1457 /* sync bCPT and fc record-keeping */
1458 if (bCPT && MASTER(cr))
1460 fcRequestCheckPoint();
1464 if (mdof_flags != 0)
1466 wallcycle_start(wcycle, ewcTRAJ);
1469 if (state->flags & (1<<estLD_RNG))
1471 get_stochd_state(upd, state);
1473 if (state->flags & (1<<estMC_RNG))
1475 get_mc_state(mcrng, state);
1481 state_global->ekinstate.bUpToDate = FALSE;
1485 update_ekinstate(&state_global->ekinstate, ekind);
1486 state_global->ekinstate.bUpToDate = TRUE;
1488 update_energyhistory(&state_global->enerhist, mdebin);
1489 if (ir->efep != efepNO || ir->bSimTemp)
1491 state_global->fep_state = state->fep_state; /* MRS: seems kludgy. The code should be
1492 structured so this isn't necessary.
1493 Note this reassignment is only necessary
1494 for single threads.*/
1495 copy_df_history(&state_global->dfhist, &df_history);
1499 write_traj(fplog, cr, outf, mdof_flags, top_global,
1500 step, t, state, state_global, f, f_global, &n_xtc, &x_xtc);
1507 if (bLastStep && step_rel == ir->nsteps &&
1508 (Flags & MD_CONFOUT) && MASTER(cr) &&
1509 !bRerunMD && !bFFscan)
1511 /* x and v have been collected in write_traj,
1512 * because a checkpoint file will always be written
1515 fprintf(stderr, "\nWriting final coordinates.\n");
1518 /* Make molecules whole only for confout writing */
1519 do_pbc_mtop(fplog, ir->ePBC, state->box, top_global, state_global->x);
1521 write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm),
1522 *top_global->name, top_global,
1523 state_global->x, state_global->v,
1524 ir->ePBC, state->box);
1527 wallcycle_stop(wcycle, ewcTRAJ);
1530 /* kludge -- virial is lost with restart for NPT control. Must restart */
1531 if (bStartingFromCpt && bVV)
1533 copy_mat(state->svir_prev, shake_vir);
1534 copy_mat(state->fvir_prev, force_vir);
1536 /* ################## END TRAJECTORY OUTPUT ################ */
1538 /* Determine the wallclock run time up till now */
1539 run_time = gmx_gettime() - (double)runtime->real;
1541 /* Check whether everything is still allright */
1542 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1543 #ifdef GMX_THREAD_MPI
1548 /* this is just make gs.sig compatible with the hack
1549 of sending signals around by MPI_Reduce with together with
1551 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1553 gs.sig[eglsSTOPCOND] = 1;
1555 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1557 gs.sig[eglsSTOPCOND] = -1;
1559 /* < 0 means stop at next step, > 0 means stop at next NS step */
1563 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1564 gmx_get_signal_name(),
1565 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1569 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1570 gmx_get_signal_name(),
1571 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1573 handled_stop_condition = (int)gmx_get_stop_condition();
1575 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1576 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1577 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1579 /* Signal to terminate the run */
1580 gs.sig[eglsSTOPCOND] = 1;
1583 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1585 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1588 if (bResetCountersHalfMaxH && MASTER(cr) &&
1589 run_time > max_hours*60.0*60.0*0.495)
1591 gs.sig[eglsRESETCOUNTERS] = 1;
1594 if (ir->nstlist == -1 && !bRerunMD)
1596 /* When bGStatEveryStep=FALSE, global_stat is only called
1597 * when we check the atom displacements, not at NS steps.
1598 * This means that also the bonded interaction count check is not
1599 * performed immediately after NS. Therefore a few MD steps could
1600 * be performed with missing interactions.
1601 * But wrong energies are never written to file,
1602 * since energies are only written after global_stat
1605 if (step >= nlh.step_nscheck)
1607 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1608 nlh.scale_tot, state->x);
1612 /* This is not necessarily true,
1613 * but step_nscheck is determined quite conservatively.
1619 /* In parallel we only have to check for checkpointing in steps
1620 * where we do global communication,
1621 * otherwise the other nodes don't know.
1623 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1626 run_time >= nchkpt*cpt_period*60.0)) &&
1627 gs.set[eglsCHKPT] == 0)
1629 gs.sig[eglsCHKPT] = 1;
1632 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1637 update_tcouple(fplog, step, ir, state, ekind, wcycle, upd, &MassQ, mdatoms);
1639 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1641 gmx_bool bIfRandomize;
1642 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1643 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1644 if (constr && bIfRandomize)
1646 update_constraints(fplog, step, &dvdl, ir, ekind, mdatoms,
1647 state, fr->bMolPBC, graph, f,
1648 &top->idef, tmp_vir, NULL,
1649 cr, nrnb, wcycle, upd, constr,
1650 bInitStep, TRUE, bCalcVir, vetanew);
1655 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1657 gmx_iterate_init(&iterate, TRUE);
1658 /* for iterations, we save these vectors, as we will be redoing the calculations */
1659 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1662 bFirstIterate = TRUE;
1663 while (bFirstIterate || iterate.bIterationActive)
1665 /* We now restore these vectors to redo the calculation with improved extended variables */
1666 if (iterate.bIterationActive)
1668 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1671 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1672 so scroll down for that logic */
1674 /* ######### START SECOND UPDATE STEP ################# */
1675 /* Box is changed in update() when we do pressure coupling,
1676 * but we should still use the old box for energy corrections and when
1677 * writing it to the energy file, so it matches the trajectory files for
1678 * the same timestep above. Make a copy in a separate array.
1680 copy_mat(state->box, lastbox);
1683 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1685 wallcycle_start(wcycle, ewcUPDATE);
1687 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1690 if (iterate.bIterationActive)
1698 /* we use a new value of scalevir to converge the iterations faster */
1699 scalevir = tracevir/trace(shake_vir);
1701 msmul(shake_vir, scalevir, shake_vir);
1702 m_add(force_vir, shake_vir, total_vir);
1703 clear_mat(shake_vir);
1705 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1706 /* We can only do Berendsen coupling after we have summed
1707 * the kinetic energy or virial. Since the happens
1708 * in global_state after update, we should only do it at
1709 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1714 update_tcouple(fplog, step, ir, state, ekind, wcycle, upd, &MassQ, mdatoms);
1715 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, wcycle,
1721 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1723 /* velocity half-step update */
1724 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1725 bUpdateDoLR, fr->f_twin, fcd,
1726 ekind, M, wcycle, upd, FALSE, etrtVELOCITY2,
1727 cr, nrnb, constr, &top->idef);
1730 /* Above, initialize just copies ekinh into ekin,
1731 * it doesn't copy position (for VV),
1732 * and entire integrator for MD.
1735 if (ir->eI == eiVVAK)
1737 copy_rvecn(state->x, cbuf, 0, state->natoms);
1739 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1741 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1742 bUpdateDoLR, fr->f_twin, fcd,
1743 ekind, M, wcycle, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1744 wallcycle_stop(wcycle, ewcUPDATE);
1746 update_constraints(fplog, step, &dvdl, ir, ekind, mdatoms, state,
1747 fr->bMolPBC, graph, f,
1748 &top->idef, shake_vir, force_vir,
1749 cr, nrnb, wcycle, upd, constr,
1750 bInitStep, FALSE, bCalcVir, state->veta);
1752 if (ir->eI == eiVVAK)
1754 /* erase F_EKIN and F_TEMP here? */
1755 /* just compute the kinetic energy at the half step to perform a trotter step */
1756 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1757 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1758 constr, NULL, FALSE, lastbox,
1759 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1760 cglo_flags | CGLO_TEMPERATURE
1762 wallcycle_start(wcycle, ewcUPDATE);
1763 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1764 /* now we know the scaling, we can compute the positions again again */
1765 copy_rvecn(cbuf, state->x, 0, state->natoms);
1767 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1769 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1770 bUpdateDoLR, fr->f_twin, fcd,
1771 ekind, M, wcycle, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1772 wallcycle_stop(wcycle, ewcUPDATE);
1774 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1775 /* are the small terms in the shake_vir here due
1776 * to numerical errors, or are they important
1777 * physically? I'm thinking they are just errors, but not completely sure.
1778 * For now, will call without actually constraining, constr=NULL*/
1779 update_constraints(fplog, step, &dvdl, ir, ekind, mdatoms,
1780 state, fr->bMolPBC, graph, f,
1781 &top->idef, tmp_vir, force_vir,
1782 cr, nrnb, wcycle, upd, NULL,
1783 bInitStep, FALSE, bCalcVir,
1786 if (!bOK && !bFFscan)
1788 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1791 if (fr->bSepDVDL && fplog && do_log)
1793 fprintf(fplog, sepdvdlformat, "Constraint dV/dl", 0.0, dvdl);
1795 enerd->term[F_DVDL_BONDED] += dvdl;
1799 /* Need to unshift here */
1800 unshift_self(graph, state->box, state->x);
1805 wallcycle_start(wcycle, ewcVSITECONSTR);
1808 shift_self(graph, state->box, state->x);
1810 construct_vsites(fplog, vsite, state->x, nrnb, ir->delta_t, state->v,
1811 top->idef.iparams, top->idef.il,
1812 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1816 unshift_self(graph, state->box, state->x);
1818 wallcycle_stop(wcycle, ewcVSITECONSTR);
1821 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1822 /* With Leap-Frog we can skip compute_globals at
1823 * non-communication steps, but we need to calculate
1824 * the kinetic energy one step before communication.
1826 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1828 if (ir->nstlist == -1 && bFirstIterate)
1830 gs.sig[eglsNABNSB] = nlh.nabnsb;
1832 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1833 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1835 bFirstIterate ? &gs : NULL,
1836 (step_rel % gs.nstms == 0) &&
1837 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1839 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1841 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1842 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1843 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1844 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1845 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1846 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1849 if (ir->nstlist == -1 && bFirstIterate)
1851 nlh.nabnsb = gs.set[eglsNABNSB];
1852 gs.set[eglsNABNSB] = 0;
1855 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1856 /* ############# END CALC EKIN AND PRESSURE ################# */
1858 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1859 the virial that should probably be addressed eventually. state->veta has better properies,
1860 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1861 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1863 if (iterate.bIterationActive &&
1864 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1865 trace(shake_vir), &tracevir))
1869 bFirstIterate = FALSE;
1872 /* only add constraint dvdl after constraints */
1873 enerd->term[F_DVDL_BONDED] += dvdl;
1874 if (!bVV || bRerunMD)
1876 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1877 sum_dhdl(enerd, state->lambda, ir->fepvals);
1879 update_box(fplog, step, ir, mdatoms, state, graph, f,
1880 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, wcycle, upd, bInitStep, FALSE);
1882 /* ################# END UPDATE STEP 2 ################# */
1883 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1885 /* The coordinates (x) were unshifted in update */
1886 if (bFFscan && (shellfc == NULL || bConverged))
1888 if (print_forcefield(fplog, enerd->term, mdatoms->homenr,
1890 &(top_global->mols), mdatoms->massT, pres))
1894 fprintf(stderr, "\n");
1900 /* We will not sum ekinh_old,
1901 * so signal that we still have to do it.
1903 bSumEkinhOld = TRUE;
1908 /* Only do GCT when the relaxation of shells (minimization) has converged,
1909 * otherwise we might be coupling to bogus energies.
1910 * In parallel we must always do this, because the other sims might
1914 /* Since this is called with the new coordinates state->x, I assume
1915 * we want the new box state->box too. / EL 20040121
1917 do_coupling(fplog, oenv, nfile, fnm, tcr, t, step, enerd->term, fr,
1919 mdatoms, &(top->idef), mu_aver,
1920 top_global->mols.nr, cr,
1921 state->box, total_vir, pres,
1922 mu_tot, state->x, f, bConverged);
1926 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1928 /* use the directly determined last velocity, not actually the averaged half steps */
1929 if (bTrotter && ir->eI == eiVV)
1931 enerd->term[F_EKIN] = last_ekin;
1933 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1937 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1941 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1943 /* Check for excessively large energies */
1947 real etot_max = 1e200;
1949 real etot_max = 1e30;
1951 if (fabs(enerd->term[F_ETOT]) > etot_max)
1953 fprintf(stderr, "Energy too large (%g), giving up\n",
1954 enerd->term[F_ETOT]);
1957 /* ######### END PREPARING EDR OUTPUT ########### */
1959 /* Time for performance */
1960 if (((step % stepout) == 0) || bLastStep)
1962 runtime_upd_proc(runtime);
1968 gmx_bool do_dr, do_or;
1970 if (fplog && do_log && bDoExpanded)
1972 /* only needed if doing expanded ensemble */
1973 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1974 &df_history, state->fep_state, ir->nstlog, step);
1976 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1980 upd_mdebin(mdebin, bDoDHDL, TRUE,
1981 t, mdatoms->tmass, enerd, state,
1982 ir->fepvals, ir->expandedvals, lastbox,
1983 shake_vir, force_vir, total_vir, pres,
1984 ekind, mu_tot, constr);
1988 upd_mdebin_step(mdebin);
1991 do_dr = do_per_step(step, ir->nstdisreout);
1992 do_or = do_per_step(step, ir->nstorireout);
1994 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
1996 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1998 if (ir->ePull != epullNO)
2000 pull_print_output(ir->pull, step, t);
2003 if (do_per_step(step, ir->nstlog))
2005 if (fflush(fplog) != 0)
2007 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
2013 /* Have to do this part after outputting the logfile and the edr file */
2014 state->fep_state = lamnew;
2015 for (i = 0; i < efptNR; i++)
2017 state_global->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
2020 /* Remaining runtime */
2021 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
2025 fprintf(stderr, "\n");
2027 print_time(stderr, runtime, step, ir, cr);
2030 /* Replica exchange */
2032 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
2033 do_per_step(step, repl_ex_nst))
2035 bExchanged = replica_exchange(fplog, cr, repl_ex,
2036 state_global, enerd,
2039 if (bExchanged && DOMAINDECOMP(cr))
2041 dd_partition_system(fplog, step, cr, TRUE, 1,
2042 state_global, top_global, ir,
2043 state, &f, mdatoms, top, fr,
2044 vsite, shellfc, constr,
2045 nrnb, wcycle, FALSE);
2051 bStartingFromCpt = FALSE;
2053 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
2054 /* With all integrators, except VV, we need to retain the pressure
2055 * at the current step for coupling at the next step.
2057 if ((state->flags & (1<<estPRES_PREV)) &&
2059 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
2061 /* Store the pressure in t_state for pressure coupling
2062 * at the next MD step.
2064 copy_mat(pres, state->pres_prev);
2067 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2069 if ( (membed != NULL) && (!bLastStep) )
2071 rescale_membed(step_rel, membed, state_global->x);
2078 /* read next frame from input trajectory */
2079 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
2084 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
2088 if (!bRerunMD || !rerun_fr.bStep)
2090 /* increase the MD step number */
2095 cycles = wallcycle_stop(wcycle, ewcSTEP);
2096 if (DOMAINDECOMP(cr) && wcycle)
2098 dd_cycles_add(cr->dd, cycles, ddCyclStep);
2101 if (bPMETuneRunning || bPMETuneTry)
2103 /* PME grid + cut-off optimization with GPUs or PME nodes */
2105 /* Count the total cycles over the last steps */
2106 cycles_pmes += cycles;
2108 /* We can only switch cut-off at NS steps */
2109 if (step % ir->nstlist == 0)
2111 /* PME grid + cut-off optimization with GPUs or PME nodes */
2114 if (DDMASTER(cr->dd))
2116 /* PME node load is too high, start tuning */
2117 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2119 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
2121 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2123 bPMETuneTry = FALSE;
2126 if (bPMETuneRunning)
2128 /* init_step might not be a multiple of nstlist,
2129 * but the first cycle is always skipped anyhow.
2132 pme_load_balance(pme_loadbal, cr,
2133 (bVerbose && MASTER(cr)) ? stderr : NULL,
2135 ir, state, cycles_pmes,
2136 fr->ic, fr->nbv, &fr->pmedata,
2139 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2140 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2141 fr->rlist = fr->ic->rlist;
2142 fr->rlistlong = fr->ic->rlistlong;
2143 fr->rcoulomb = fr->ic->rcoulomb;
2144 fr->rvdw = fr->ic->rvdw;
2150 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2151 gs.set[eglsRESETCOUNTERS] != 0)
2153 /* Reset all the counters related to performance over the run */
2154 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, runtime,
2155 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2156 wcycle_set_reset_counters(wcycle, -1);
2157 if (!(cr->duty & DUTY_PME))
2159 /* Tell our PME node to reset its counters */
2160 gmx_pme_send_resetcounters(cr, step);
2162 /* Correct max_hours for the elapsed time */
2163 max_hours -= run_time/(60.0*60.0);
2164 bResetCountersHalfMaxH = FALSE;
2165 gs.set[eglsRESETCOUNTERS] = 0;
2169 /* End of main MD loop */
2173 runtime_end(runtime);
2175 if (bRerunMD && MASTER(cr))
2180 if (!(cr->duty & DUTY_PME))
2182 /* Tell the PME only node to finish */
2183 gmx_pme_send_finish(cr);
2188 if (ir->nstcalcenergy > 0 && !bRerunMD)
2190 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
2191 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
2199 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2201 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)));
2202 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
2205 if (pme_loadbal != NULL)
2207 pme_loadbal_done(pme_loadbal, fplog);
2210 if (shellfc && fplog)
2212 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
2213 (nconverged*100.0)/step_rel);
2214 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
2218 if (repl_ex_nst > 0 && MASTER(cr))
2220 print_replica_exchange_statistics(fplog, repl_ex);
2223 runtime->nsteps_done = step_rel;