3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
55 #include "md_support.h"
56 #include "md_logging.h"
71 #include "domdec_network.h"
77 #include "compute_io.h"
79 #include "checkpoint.h"
80 #include "mtop_util.h"
81 #include "sighandler.h"
84 #include "pme_loadbal.h"
87 #include "types/nlistheuristics.h"
88 #include "types/iteratedconstraints.h"
89 #include "nbnxn_cuda_data_mgmt.h"
91 #include "gromacs/utility/gmxmpi.h"
97 static void reset_all_counters(FILE *fplog, t_commrec *cr,
99 gmx_large_int_t *step_rel, t_inputrec *ir,
100 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
101 gmx_runtime_t *runtime,
102 nbnxn_cuda_ptr_t cu_nbv)
104 char sbuf[STEPSTRSIZE];
106 /* Reset all the counters related to performance over the run */
107 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
108 gmx_step_str(step, sbuf));
112 nbnxn_cuda_reset_timings(cu_nbv);
115 wallcycle_stop(wcycle, ewcRUN);
116 wallcycle_reset_all(wcycle);
117 if (DOMAINDECOMP(cr))
119 reset_dd_statistics_counters(cr->dd);
122 ir->init_step += *step_rel;
123 ir->nsteps -= *step_rel;
125 wallcycle_start(wcycle, ewcRUN);
126 runtime_start(runtime);
127 print_date_and_time(fplog, cr->nodeid, "Restarted time", runtime);
130 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
131 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
133 gmx_vsite_t *vsite, gmx_constr_t constr,
134 int stepout, t_inputrec *ir,
135 gmx_mtop_t *top_global,
137 t_state *state_global,
139 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
140 gmx_edsam_t ed, t_forcerec *fr,
141 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
142 real cpt_period, real max_hours,
143 const char gmx_unused *deviceOptions,
145 gmx_runtime_t *runtime)
148 gmx_large_int_t step, step_rel;
150 double t, t0, lam0[efptNR];
151 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
152 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
153 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
154 bBornRadii, bStartingFromCpt;
155 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
156 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
157 bForceUpdate = FALSE, bCPT;
158 gmx_bool bMasterState;
159 int force_flags, cglo_flags;
160 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
165 t_state *bufstate = NULL;
166 matrix *scale_tot, pcoupl_mu, M, ebox;
169 gmx_repl_ex_t repl_ex = NULL;
172 t_mdebin *mdebin = NULL;
173 df_history_t df_history;
174 t_state *state = NULL;
175 rvec *f_global = NULL;
176 gmx_enerdata_t *enerd;
178 gmx_global_stat_t gstat;
179 gmx_update_t upd = NULL;
180 t_graph *graph = NULL;
182 gmx_rng_t mcrng = NULL;
183 gmx_groups_t *groups;
184 gmx_ekindata_t *ekind, *ekind_save;
185 gmx_shellfc_t shellfc;
186 int count, nconverged = 0;
189 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
191 gmx_bool bResetCountersHalfMaxH = FALSE;
192 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
193 gmx_bool bUpdateDoLR;
198 real veta_save, scalevir, tracevir;
204 real saved_conserved_quantity = 0;
209 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
210 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
211 gmx_iterate_t iterate;
212 gmx_large_int_t multisim_nsteps = -1; /* number of steps to do before first multisim
213 simulation stops. If equal to zero, don't
214 communicate any more between multisims.*/
215 /* PME load balancing data for GPU kernels */
216 pme_load_balancing_t pme_loadbal = NULL;
218 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
221 /* Temporary addition for FAHCORE checkpointing */
225 /* Check for special mdrun options */
226 bRerunMD = (Flags & MD_RERUN);
227 bAppend = (Flags & MD_APPENDFILES);
228 if (Flags & MD_RESETCOUNTERSHALFWAY)
232 /* Signal to reset the counters half the simulation steps. */
233 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
235 /* Signal to reset the counters halfway the simulation time. */
236 bResetCountersHalfMaxH = (max_hours > 0);
239 /* md-vv uses averaged full step velocities for T-control
240 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
241 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
243 if (bVV) /* to store the initial velocities while computing virial */
245 snew(cbuf, top_global->natoms);
247 /* all the iteratative cases - only if there are constraints */
248 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
249 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
250 false in this step. The correct value, true or false,
251 is set at each step, as it depends on the frequency of temperature
252 and pressure control.*/
253 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
257 /* Since we don't know if the frames read are related in any way,
258 * rebuild the neighborlist at every step.
261 ir->nstcalcenergy = 1;
265 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
267 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
268 bGStatEveryStep = (nstglobalcomm == 1);
270 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
273 "To reduce the energy communication with nstlist = -1\n"
274 "the neighbor list validity should not be checked at every step,\n"
275 "this means that exact integration is not guaranteed.\n"
276 "The neighbor list validity is checked after:\n"
277 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
278 "In most cases this will result in exact integration.\n"
279 "This reduces the energy communication by a factor of 2 to 3.\n"
280 "If you want less energy communication, set nstlist > 3.\n\n");
287 groups = &top_global->groups;
290 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
291 &(state_global->fep_state), lam0,
292 nrnb, top_global, &upd,
293 nfile, fnm, &outf, &mdebin,
294 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
296 clear_mat(total_vir);
298 /* Energy terms and groups */
300 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
302 if (DOMAINDECOMP(cr))
308 snew(f, top_global->natoms);
311 /* lambda Monte carlo random number generator */
314 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
316 /* copy the state into df_history */
317 copy_df_history(&df_history, &state_global->dfhist);
319 /* Kinetic energy data */
321 init_ekindata(fplog, top_global, &(ir->opts), ekind);
322 /* needed for iteration of constraints */
324 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
325 /* Copy the cos acceleration to the groups struct */
326 ekind->cosacc.cos_accel = ir->cos_accel;
328 gstat = global_stat_init(ir);
331 /* Check for polarizable models and flexible constraints */
332 shellfc = init_shell_flexcon(fplog,
333 top_global, n_flexible_constraints(constr),
334 (ir->bContinuation ||
335 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
336 NULL : state_global->x);
340 #ifdef GMX_THREAD_MPI
341 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
343 set_deform_reference_box(upd,
344 deform_init_init_step_tpx,
345 deform_init_box_tpx);
346 #ifdef GMX_THREAD_MPI
347 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
352 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
353 if ((io > 2000) && MASTER(cr))
356 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
361 if (DOMAINDECOMP(cr))
363 top = dd_init_local_top(top_global);
366 dd_init_local_state(cr->dd, state_global, state);
368 if (DDMASTER(cr->dd) && ir->nstfout)
370 snew(f_global, state_global->natoms);
377 /* Initialize the particle decomposition and split the topology */
378 top = split_system(fplog, top_global, ir, cr);
380 pd_cg_range(cr, &fr->cg0, &fr->hcg);
381 pd_at_range(cr, &a0, &a1);
385 top = gmx_mtop_generate_local_top(top_global, ir);
388 a1 = top_global->natoms;
391 forcerec_set_excl_load(fr, top, cr);
393 state = partdec_init_local_state(cr, state_global);
396 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
400 set_vsite_top(vsite, top, mdatoms, cr);
403 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
405 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
410 make_local_shells(cr, mdatoms, shellfc);
413 setup_bonded_threading(fr, &top->idef);
415 if (ir->pull && PAR(cr))
417 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
421 if (DOMAINDECOMP(cr))
423 /* Distribute the charge groups over the nodes from the master node */
424 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
425 state_global, top_global, ir,
426 state, &f, mdatoms, top, fr,
427 vsite, shellfc, constr,
428 nrnb, wcycle, FALSE);
432 update_mdatoms(mdatoms, state->lambda[efptMASS]);
434 if (opt2bSet("-cpi", nfile, fnm))
436 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
440 bStateFromCP = FALSE;
447 /* Update mdebin with energy history if appending to output files */
448 if (Flags & MD_APPENDFILES)
450 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
454 /* We might have read an energy history from checkpoint,
455 * free the allocated memory and reset the counts.
457 done_energyhistory(&state_global->enerhist);
458 init_energyhistory(&state_global->enerhist);
461 /* Set the initial energy history in state by updating once */
462 update_energyhistory(&state_global->enerhist, mdebin);
465 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
467 /* Set the random state if we read a checkpoint file */
468 set_stochd_state(upd, state);
471 if (state->flags & (1<<estMC_RNG))
473 set_mc_state(mcrng, state);
476 /* Initialize constraints */
479 if (!DOMAINDECOMP(cr))
481 set_constraints(constr, top, ir, mdatoms, cr);
487 /* We need to be sure replica exchange can only occur
488 * when the energies are current */
489 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
490 "repl_ex_nst", &repl_ex_nst);
491 /* This check needs to happen before inter-simulation
492 * signals are initialized, too */
494 if (repl_ex_nst > 0 && MASTER(cr))
496 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
497 repl_ex_nst, repl_ex_nex, repl_ex_seed);
500 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
501 * With perturbed charges with soft-core we should not change the cut-off.
503 if ((Flags & MD_TUNEPME) &&
504 EEL_PME(fr->eeltype) &&
505 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
506 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
509 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
511 if (cr->duty & DUTY_PME)
513 /* Start tuning right away, as we can't measure the load */
514 bPMETuneRunning = TRUE;
518 /* Separate PME nodes, we can measure the PP/PME load balance */
523 if (!ir->bContinuation && !bRerunMD)
525 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
527 /* Set the velocities of frozen particles to zero */
528 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
530 for (m = 0; m < DIM; m++)
532 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
542 /* Constrain the initial coordinates and velocities */
543 do_constrain_first(fplog, constr, ir, mdatoms, state,
548 /* Construct the virtual sites for the initial configuration */
549 construct_vsites(vsite, state->x, ir->delta_t, NULL,
550 top->idef.iparams, top->idef.il,
551 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
557 /* set free energy calculation frequency as the minimum
558 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
559 nstfep = ir->fepvals->nstdhdl;
562 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
566 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
569 /* I'm assuming we need global communication the first time! MRS */
570 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
571 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
572 | (bVV ? CGLO_PRESSURE : 0)
573 | (bVV ? CGLO_CONSTRAINT : 0)
574 | (bRerunMD ? CGLO_RERUNMD : 0)
575 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
577 bSumEkinhOld = FALSE;
578 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
579 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
580 constr, NULL, FALSE, state->box,
581 top_global, &bSumEkinhOld, cglo_flags);
582 if (ir->eI == eiVVAK)
584 /* a second call to get the half step temperature initialized as well */
585 /* we do the same call as above, but turn the pressure off -- internally to
586 compute_globals, this is recognized as a velocity verlet half-step
587 kinetic energy calculation. This minimized excess variables, but
588 perhaps loses some logic?*/
590 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
591 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
592 constr, NULL, FALSE, state->box,
593 top_global, &bSumEkinhOld,
594 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
597 /* Calculate the initial half step temperature, and save the ekinh_old */
598 if (!(Flags & MD_STARTFROMCPT))
600 for (i = 0; (i < ir->opts.ngtc); i++)
602 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
607 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
608 and there is no previous step */
611 /* if using an iterative algorithm, we need to create a working directory for the state. */
614 bufstate = init_bufstate(state);
617 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
618 temperature control */
619 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
623 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
626 "RMS relative constraint deviation after constraining: %.2e\n",
627 constr_rmsd(constr, FALSE));
629 if (EI_STATE_VELOCITY(ir->eI))
631 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
635 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
636 " input trajectory '%s'\n\n",
637 *(top_global->name), opt2fn("-rerun", nfile, fnm));
640 fprintf(stderr, "Calculated time to finish depends on nsteps from "
641 "run input file,\nwhich may not correspond to the time "
642 "needed to process input trajectory.\n\n");
648 fprintf(stderr, "starting mdrun '%s'\n",
649 *(top_global->name));
652 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
656 sprintf(tbuf, "%s", "infinite");
658 if (ir->init_step > 0)
660 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
661 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
662 gmx_step_str(ir->init_step, sbuf2),
663 ir->init_step*ir->delta_t);
667 fprintf(stderr, "%s steps, %s ps.\n",
668 gmx_step_str(ir->nsteps, sbuf), tbuf);
671 fprintf(fplog, "\n");
674 /* Set and write start time */
675 runtime_start(runtime);
676 print_date_and_time(fplog, cr->nodeid, "Started mdrun", runtime);
677 wallcycle_start(wcycle, ewcRUN);
680 fprintf(fplog, "\n");
683 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
685 chkpt_ret = fcCheckPointParallel( cr->nodeid,
689 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
694 /***********************************************************
698 ************************************************************/
700 /* if rerunMD then read coordinates and velocities from input trajectory */
703 if (getenv("GMX_FORCE_UPDATE"))
711 bNotLastFrame = read_first_frame(oenv, &status,
712 opt2fn("-rerun", nfile, fnm),
713 &rerun_fr, TRX_NEED_X | TRX_READ_V);
714 if (rerun_fr.natoms != top_global->natoms)
717 "Number of atoms in trajectory (%d) does not match the "
718 "run input file (%d)\n",
719 rerun_fr.natoms, top_global->natoms);
721 if (ir->ePBC != epbcNONE)
725 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);
727 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
729 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
736 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
739 if (ir->ePBC != epbcNONE)
741 /* Set the shift vectors.
742 * Necessary here when have a static box different from the tpr box.
744 calc_shifts(rerun_fr.box, fr->shift_vec);
748 /* loop over MD steps or if rerunMD to end of input trajectory */
750 /* Skip the first Nose-Hoover integration when we get the state from tpx */
751 bStateFromTPX = !bStateFromCP;
752 bInitStep = bFirstStep && (bStateFromTPX || bVV);
753 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
755 bSumEkinhOld = FALSE;
758 init_global_signals(&gs, cr, ir, repl_ex_nst);
760 step = ir->init_step;
763 if (ir->nstlist == -1)
765 init_nlistheuristics(&nlh, bGStatEveryStep, step);
768 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
770 /* check how many steps are left in other sims */
771 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
775 /* and stop now if we should */
776 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
777 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
778 while (!bLastStep || (bRerunMD && bNotLastFrame))
781 wallcycle_start(wcycle, ewcSTEP);
787 step = rerun_fr.step;
788 step_rel = step - ir->init_step;
801 bLastStep = (step_rel == ir->nsteps);
802 t = t0 + step*ir->delta_t;
805 if (ir->efep != efepNO || ir->bSimTemp)
807 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
808 requiring different logic. */
810 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
811 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
812 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
813 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
818 update_annealing_target_temp(&(ir->opts), t);
823 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
825 for (i = 0; i < state_global->natoms; i++)
827 copy_rvec(rerun_fr.x[i], state_global->x[i]);
831 for (i = 0; i < state_global->natoms; i++)
833 copy_rvec(rerun_fr.v[i], state_global->v[i]);
838 for (i = 0; i < state_global->natoms; i++)
840 clear_rvec(state_global->v[i]);
844 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
845 " Ekin, temperature and pressure are incorrect,\n"
846 " the virial will be incorrect when constraints are present.\n"
848 bRerunWarnNoV = FALSE;
852 copy_mat(rerun_fr.box, state_global->box);
853 copy_mat(state_global->box, state->box);
855 if (vsite && (Flags & MD_RERUN_VSITE))
857 if (DOMAINDECOMP(cr))
859 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
863 /* Following is necessary because the graph may get out of sync
864 * with the coordinates if we only have every N'th coordinate set
866 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
867 shift_self(graph, state->box, state->x);
869 construct_vsites(vsite, state->x, ir->delta_t, state->v,
870 top->idef.iparams, top->idef.il,
871 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
874 unshift_self(graph, state->box, state->x);
879 /* Stop Center of Mass motion */
880 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
884 /* for rerun MD always do Neighbour Searching */
885 bNS = (bFirstStep || ir->nstlist != 0);
890 /* Determine whether or not to do Neighbour Searching and LR */
891 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
893 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
894 (ir->nstlist == -1 && nlh.nabnsb > 0));
896 if (bNS && ir->nstlist == -1)
898 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
902 /* check whether we should stop because another simulation has
906 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
907 (multisim_nsteps != ir->nsteps) )
914 "Stopping simulation %d because another one has finished\n",
918 gs.sig[eglsCHKPT] = 1;
923 /* < 0 means stop at next step, > 0 means stop at next NS step */
924 if ( (gs.set[eglsSTOPCOND] < 0) ||
925 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
930 /* Determine whether or not to update the Born radii if doing GB */
931 bBornRadii = bFirstStep;
932 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
937 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
938 do_verbose = bVerbose &&
939 (step % stepout == 0 || bFirstStep || bLastStep);
941 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
949 bMasterState = FALSE;
950 /* Correct the new box if it is too skewed */
951 if (DYNAMIC_BOX(*ir))
953 if (correct_box(fplog, step, state->box, graph))
958 if (DOMAINDECOMP(cr) && bMasterState)
960 dd_collect_state(cr->dd, state, state_global);
964 if (DOMAINDECOMP(cr))
966 /* Repartition the domain decomposition */
967 wallcycle_start(wcycle, ewcDOMDEC);
968 dd_partition_system(fplog, step, cr,
969 bMasterState, nstglobalcomm,
970 state_global, top_global, ir,
971 state, &f, mdatoms, top, fr,
972 vsite, shellfc, constr,
974 do_verbose && !bPMETuneRunning);
975 wallcycle_stop(wcycle, ewcDOMDEC);
976 /* If using an iterative integrator, reallocate space to match the decomposition */
980 if (MASTER(cr) && do_log)
982 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
985 if (ir->efep != efepNO)
987 update_mdatoms(mdatoms, state->lambda[efptMASS]);
990 if ((bRerunMD && rerun_fr.bV) || bExchanged)
993 /* We need the kinetic energy at minus the half step for determining
994 * the full step kinetic energy and possibly for T-coupling.*/
995 /* This may not be quite working correctly yet . . . . */
996 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
997 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
998 constr, NULL, FALSE, state->box,
999 top_global, &bSumEkinhOld,
1000 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1002 clear_mat(force_vir);
1004 /* We write a checkpoint at this MD step when:
1005 * either at an NS step when we signalled through gs,
1006 * or at the last step (but not when we do not want confout),
1007 * but never at the first step or with rerun.
1009 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1010 (bLastStep && (Flags & MD_CONFOUT))) &&
1011 step > ir->init_step && !bRerunMD);
1014 gs.set[eglsCHKPT] = 0;
1017 /* Determine the energy and pressure:
1018 * at nstcalcenergy steps and at energy output steps (set below).
1020 if (EI_VV(ir->eI) && (!bInitStep))
1022 /* for vv, the first half of the integration actually corresponds
1023 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1024 but the virial needs to be calculated on both the current step and the 'next' step. Future
1025 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1027 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1028 bCalcVir = bCalcEner ||
1029 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1033 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1034 bCalcVir = bCalcEner ||
1035 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1038 /* Do we need global communication ? */
1039 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1040 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1041 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1043 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1045 if (do_ene || do_log)
1052 /* these CGLO_ options remain the same throughout the iteration */
1053 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1054 (bGStat ? CGLO_GSTAT : 0)
1057 force_flags = (GMX_FORCE_STATECHANGED |
1058 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1059 GMX_FORCE_ALLFORCES |
1061 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1062 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1063 (bDoFEP ? GMX_FORCE_DHDL : 0)
1068 if (do_per_step(step, ir->nstcalclr))
1070 force_flags |= GMX_FORCE_DO_LR;
1076 /* Now is the time to relax the shells */
1077 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1078 ir, bNS, force_flags,
1081 state, f, force_vir, mdatoms,
1082 nrnb, wcycle, graph, groups,
1083 shellfc, fr, bBornRadii, t, mu_tot,
1095 /* The coordinates (x) are shifted (to get whole molecules)
1097 * This is parallellized as well, and does communication too.
1098 * Check comments in sim_util.c
1100 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1101 state->box, state->x, &state->hist,
1102 f, force_vir, mdatoms, enerd, fcd,
1103 state->lambda, graph,
1104 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1105 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1108 if (bVV && !bStartingFromCpt && !bRerunMD)
1109 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1111 if (ir->eI == eiVV && bInitStep)
1113 /* if using velocity verlet with full time step Ekin,
1114 * take the first half step only to compute the
1115 * virial for the first step. From there,
1116 * revert back to the initial coordinates
1117 * so that the input is actually the initial step.
1119 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1123 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1124 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1127 /* If we are using twin-range interactions where the long-range component
1128 * is only evaluated every nstcalclr>1 steps, we should do a special update
1129 * step to combine the long-range forces on these steps.
1130 * For nstcalclr=1 this is not done, since the forces would have been added
1131 * directly to the short-range forces already.
1133 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1135 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1136 f, bUpdateDoLR, fr->f_twin, fcd,
1137 ekind, M, upd, bInitStep, etrtVELOCITY1,
1138 cr, nrnb, constr, &top->idef);
1140 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1142 gmx_iterate_init(&iterate, TRUE);
1144 /* for iterations, we save these vectors, as we will be self-consistently iterating
1147 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1149 /* save the state */
1150 if (iterate.bIterationActive)
1152 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1155 bFirstIterate = TRUE;
1156 while (bFirstIterate || iterate.bIterationActive)
1158 if (iterate.bIterationActive)
1160 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1161 if (bFirstIterate && bTrotter)
1163 /* The first time through, we need a decent first estimate
1164 of veta(t+dt) to compute the constraints. Do
1165 this by computing the box volume part of the
1166 trotter integration at this time. Nothing else
1167 should be changed by this routine here. If
1168 !(first time), we start with the previous value
1171 veta_save = state->veta;
1172 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1173 vetanew = state->veta;
1174 state->veta = veta_save;
1179 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1181 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1182 state, fr->bMolPBC, graph, f,
1183 &top->idef, shake_vir,
1184 cr, nrnb, wcycle, upd, constr,
1185 TRUE, bCalcVir, vetanew);
1189 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1195 /* Need to unshift here if a do_force has been
1196 called in the previous step */
1197 unshift_self(graph, state->box, state->x);
1200 /* if VV, compute the pressure and constraints */
1201 /* For VV2, we strictly only need this if using pressure
1202 * control, but we really would like to have accurate pressures
1204 * Think about ways around this in the future?
1205 * For now, keep this choice in comments.
1207 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1208 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1210 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1211 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1213 bSumEkinhOld = TRUE;
1215 /* for vv, the first half of the integration actually corresponds to the previous step.
1216 So we need information from the last step in the first half of the integration */
1217 if (bGStat || do_per_step(step-1, nstglobalcomm))
1219 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1220 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1221 constr, NULL, FALSE, state->box,
1222 top_global, &bSumEkinhOld,
1225 | (bTemp ? CGLO_TEMPERATURE : 0)
1226 | (bPres ? CGLO_PRESSURE : 0)
1227 | (bPres ? CGLO_CONSTRAINT : 0)
1228 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1229 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1232 /* explanation of above:
1233 a) We compute Ekin at the full time step
1234 if 1) we are using the AveVel Ekin, and it's not the
1235 initial step, or 2) if we are using AveEkin, but need the full
1236 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1237 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1238 EkinAveVel because it's needed for the pressure */
1240 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1245 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1246 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1253 /* We need the kinetic energy at minus the half step for determining
1254 * the full step kinetic energy and possibly for T-coupling.*/
1255 /* This may not be quite working correctly yet . . . . */
1256 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1257 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1258 constr, NULL, FALSE, state->box,
1259 top_global, &bSumEkinhOld,
1260 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1265 if (iterate.bIterationActive &&
1266 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1267 state->veta, &vetanew))
1271 bFirstIterate = FALSE;
1274 if (bTrotter && !bInitStep)
1276 copy_mat(shake_vir, state->svir_prev);
1277 copy_mat(force_vir, state->fvir_prev);
1278 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1280 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1281 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1282 enerd->term[F_EKIN] = trace(ekind->ekin);
1285 /* if it's the initial step, we performed this first step just to get the constraint virial */
1286 if (bInitStep && ir->eI == eiVV)
1288 copy_rvecn(cbuf, state->v, 0, state->natoms);
1292 /* MRS -- now done iterating -- compute the conserved quantity */
1295 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1298 last_ekin = enerd->term[F_EKIN];
1300 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1302 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1304 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1307 sum_dhdl(enerd, state->lambda, ir->fepvals);
1311 /* ######## END FIRST UPDATE STEP ############## */
1312 /* ######## If doing VV, we now have v(dt) ###### */
1315 /* perform extended ensemble sampling in lambda - we don't
1316 actually move to the new state before outputting
1317 statistics, but if performing simulated tempering, we
1318 do update the velocities and the tau_t. */
1320 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, &df_history, step, mcrng, state->v, mdatoms);
1323 /* Now we have the energies and forces corresponding to the
1324 * coordinates at time t. We must output all of this before
1327 do_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1328 ir, state, state_global, top_global, fr, upd,
1329 outf, mdebin, ekind, df_history, f, f_global,
1330 wcycle, mcrng, &nchkpt,
1331 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1334 /* kludge -- virial is lost with restart for NPT control. Must restart */
1335 if (bStartingFromCpt && bVV)
1337 copy_mat(state->svir_prev, shake_vir);
1338 copy_mat(state->fvir_prev, force_vir);
1341 /* Determine the wallclock run time up till now */
1342 run_time = gmx_gettime() - (double)runtime->real;
1344 /* Check whether everything is still allright */
1345 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1346 #ifdef GMX_THREAD_MPI
1351 /* this is just make gs.sig compatible with the hack
1352 of sending signals around by MPI_Reduce with together with
1354 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1356 gs.sig[eglsSTOPCOND] = 1;
1358 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1360 gs.sig[eglsSTOPCOND] = -1;
1362 /* < 0 means stop at next step, > 0 means stop at next NS step */
1366 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1367 gmx_get_signal_name(),
1368 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1372 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1373 gmx_get_signal_name(),
1374 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1376 handled_stop_condition = (int)gmx_get_stop_condition();
1378 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1379 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1380 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1382 /* Signal to terminate the run */
1383 gs.sig[eglsSTOPCOND] = 1;
1386 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1388 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1391 if (bResetCountersHalfMaxH && MASTER(cr) &&
1392 run_time > max_hours*60.0*60.0*0.495)
1394 gs.sig[eglsRESETCOUNTERS] = 1;
1397 if (ir->nstlist == -1 && !bRerunMD)
1399 /* When bGStatEveryStep=FALSE, global_stat is only called
1400 * when we check the atom displacements, not at NS steps.
1401 * This means that also the bonded interaction count check is not
1402 * performed immediately after NS. Therefore a few MD steps could
1403 * be performed with missing interactions.
1404 * But wrong energies are never written to file,
1405 * since energies are only written after global_stat
1408 if (step >= nlh.step_nscheck)
1410 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1411 nlh.scale_tot, state->x);
1415 /* This is not necessarily true,
1416 * but step_nscheck is determined quite conservatively.
1422 /* In parallel we only have to check for checkpointing in steps
1423 * where we do global communication,
1424 * otherwise the other nodes don't know.
1426 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1429 run_time >= nchkpt*cpt_period*60.0)) &&
1430 gs.set[eglsCHKPT] == 0)
1432 gs.sig[eglsCHKPT] = 1;
1435 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1440 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1442 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1444 gmx_bool bIfRandomize;
1445 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1446 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1447 if (constr && bIfRandomize)
1449 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1450 state, fr->bMolPBC, graph, f,
1451 &top->idef, tmp_vir,
1452 cr, nrnb, wcycle, upd, constr,
1453 TRUE, bCalcVir, vetanew);
1458 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1460 gmx_iterate_init(&iterate, TRUE);
1461 /* for iterations, we save these vectors, as we will be redoing the calculations */
1462 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1465 bFirstIterate = TRUE;
1466 while (bFirstIterate || iterate.bIterationActive)
1468 /* We now restore these vectors to redo the calculation with improved extended variables */
1469 if (iterate.bIterationActive)
1471 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1474 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1475 so scroll down for that logic */
1477 /* ######### START SECOND UPDATE STEP ################# */
1478 /* Box is changed in update() when we do pressure coupling,
1479 * but we should still use the old box for energy corrections and when
1480 * writing it to the energy file, so it matches the trajectory files for
1481 * the same timestep above. Make a copy in a separate array.
1483 copy_mat(state->box, lastbox);
1488 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1490 wallcycle_start(wcycle, ewcUPDATE);
1491 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1494 if (iterate.bIterationActive)
1502 /* we use a new value of scalevir to converge the iterations faster */
1503 scalevir = tracevir/trace(shake_vir);
1505 msmul(shake_vir, scalevir, shake_vir);
1506 m_add(force_vir, shake_vir, total_vir);
1507 clear_mat(shake_vir);
1509 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1510 /* We can only do Berendsen coupling after we have summed
1511 * the kinetic energy or virial. Since the happens
1512 * in global_state after update, we should only do it at
1513 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1518 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1519 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1524 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1526 /* velocity half-step update */
1527 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1528 bUpdateDoLR, fr->f_twin, fcd,
1529 ekind, M, upd, FALSE, etrtVELOCITY2,
1530 cr, nrnb, constr, &top->idef);
1533 /* Above, initialize just copies ekinh into ekin,
1534 * it doesn't copy position (for VV),
1535 * and entire integrator for MD.
1538 if (ir->eI == eiVVAK)
1540 copy_rvecn(state->x, cbuf, 0, state->natoms);
1542 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1544 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1545 bUpdateDoLR, fr->f_twin, fcd,
1546 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1547 wallcycle_stop(wcycle, ewcUPDATE);
1549 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1550 fr->bMolPBC, graph, f,
1551 &top->idef, shake_vir,
1552 cr, nrnb, wcycle, upd, constr,
1553 FALSE, bCalcVir, state->veta);
1555 if (ir->eI == eiVVAK)
1557 /* erase F_EKIN and F_TEMP here? */
1558 /* just compute the kinetic energy at the half step to perform a trotter step */
1559 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1560 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1561 constr, NULL, FALSE, lastbox,
1562 top_global, &bSumEkinhOld,
1563 cglo_flags | CGLO_TEMPERATURE
1565 wallcycle_start(wcycle, ewcUPDATE);
1566 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1567 /* now we know the scaling, we can compute the positions again again */
1568 copy_rvecn(cbuf, state->x, 0, state->natoms);
1570 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1572 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1573 bUpdateDoLR, fr->f_twin, fcd,
1574 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1575 wallcycle_stop(wcycle, ewcUPDATE);
1577 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1578 /* are the small terms in the shake_vir here due
1579 * to numerical errors, or are they important
1580 * physically? I'm thinking they are just errors, but not completely sure.
1581 * For now, will call without actually constraining, constr=NULL*/
1582 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1583 state, fr->bMolPBC, graph, f,
1584 &top->idef, tmp_vir,
1585 cr, nrnb, wcycle, upd, NULL,
1591 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1594 if (fr->bSepDVDL && fplog && do_log)
1596 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1600 /* this factor or 2 correction is necessary
1601 because half of the constraint force is removed
1602 in the vv step, so we have to double it. See
1603 the Redmine issue #1255. It is not yet clear
1604 if the factor of 2 is exact, or just a very
1605 good approximation, and this will be
1606 investigated. The next step is to see if this
1607 can be done adding a dhdl contribution from the
1608 rattle step, but this is somewhat more
1609 complicated with the current code. Will be
1610 investigated, hopefully for 4.6.3. However,
1611 this current solution is much better than
1612 having it completely wrong.
1614 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1618 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1623 /* Need to unshift here */
1624 unshift_self(graph, state->box, state->x);
1629 wallcycle_start(wcycle, ewcVSITECONSTR);
1632 shift_self(graph, state->box, state->x);
1634 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1635 top->idef.iparams, top->idef.il,
1636 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1640 unshift_self(graph, state->box, state->x);
1642 wallcycle_stop(wcycle, ewcVSITECONSTR);
1645 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1646 /* With Leap-Frog we can skip compute_globals at
1647 * non-communication steps, but we need to calculate
1648 * the kinetic energy one step before communication.
1650 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1652 if (ir->nstlist == -1 && bFirstIterate)
1654 gs.sig[eglsNABNSB] = nlh.nabnsb;
1656 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1657 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1659 bFirstIterate ? &gs : NULL,
1660 (step_rel % gs.nstms == 0) &&
1661 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1663 top_global, &bSumEkinhOld,
1665 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1666 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1667 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1668 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1669 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1670 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1673 if (ir->nstlist == -1 && bFirstIterate)
1675 nlh.nabnsb = gs.set[eglsNABNSB];
1676 gs.set[eglsNABNSB] = 0;
1679 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1680 /* ############# END CALC EKIN AND PRESSURE ################# */
1682 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1683 the virial that should probably be addressed eventually. state->veta has better properies,
1684 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1685 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1687 if (iterate.bIterationActive &&
1688 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1689 trace(shake_vir), &tracevir))
1693 bFirstIterate = FALSE;
1696 if (!bVV || bRerunMD)
1698 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1699 sum_dhdl(enerd, state->lambda, ir->fepvals);
1701 update_box(fplog, step, ir, mdatoms, state, f,
1702 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1704 /* ################# END UPDATE STEP 2 ################# */
1705 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1707 /* The coordinates (x) were unshifted in update */
1710 /* We will not sum ekinh_old,
1711 * so signal that we still have to do it.
1713 bSumEkinhOld = TRUE;
1716 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1718 /* use the directly determined last velocity, not actually the averaged half steps */
1719 if (bTrotter && ir->eI == eiVV)
1721 enerd->term[F_EKIN] = last_ekin;
1723 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1727 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1731 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1733 /* ######### END PREPARING EDR OUTPUT ########### */
1735 /* Time for performance */
1736 if (((step % stepout) == 0) || bLastStep)
1738 runtime_upd_proc(runtime);
1744 gmx_bool do_dr, do_or;
1746 if (fplog && do_log && bDoExpanded)
1748 /* only needed if doing expanded ensemble */
1749 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1750 &df_history, state->fep_state, ir->nstlog, step);
1752 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1756 upd_mdebin(mdebin, bDoDHDL, TRUE,
1757 t, mdatoms->tmass, enerd, state,
1758 ir->fepvals, ir->expandedvals, lastbox,
1759 shake_vir, force_vir, total_vir, pres,
1760 ekind, mu_tot, constr);
1764 upd_mdebin_step(mdebin);
1767 do_dr = do_per_step(step, ir->nstdisreout);
1768 do_or = do_per_step(step, ir->nstorireout);
1770 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
1772 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1774 if (ir->ePull != epullNO)
1776 pull_print_output(ir->pull, step, t);
1779 if (do_per_step(step, ir->nstlog))
1781 if (fflush(fplog) != 0)
1783 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1789 /* Have to do this part after outputting the logfile and the edr file */
1790 state->fep_state = lamnew;
1791 for (i = 0; i < efptNR; i++)
1793 state_global->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
1796 /* Remaining runtime */
1797 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1801 fprintf(stderr, "\n");
1803 print_time(stderr, runtime, step, ir, cr);
1806 /* Replica exchange */
1808 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1809 do_per_step(step, repl_ex_nst))
1811 bExchanged = replica_exchange(fplog, cr, repl_ex,
1812 state_global, enerd,
1815 if (bExchanged && DOMAINDECOMP(cr))
1817 dd_partition_system(fplog, step, cr, TRUE, 1,
1818 state_global, top_global, ir,
1819 state, &f, mdatoms, top, fr,
1820 vsite, shellfc, constr,
1821 nrnb, wcycle, FALSE);
1827 bStartingFromCpt = FALSE;
1829 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1830 /* With all integrators, except VV, we need to retain the pressure
1831 * at the current step for coupling at the next step.
1833 if ((state->flags & (1<<estPRES_PREV)) &&
1835 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1837 /* Store the pressure in t_state for pressure coupling
1838 * at the next MD step.
1840 copy_mat(pres, state->pres_prev);
1843 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1845 if ( (membed != NULL) && (!bLastStep) )
1847 rescale_membed(step_rel, membed, state_global->x);
1854 /* read next frame from input trajectory */
1855 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1860 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1864 if (!bRerunMD || !rerun_fr.bStep)
1866 /* increase the MD step number */
1871 cycles = wallcycle_stop(wcycle, ewcSTEP);
1872 if (DOMAINDECOMP(cr) && wcycle)
1874 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1877 if (bPMETuneRunning || bPMETuneTry)
1879 /* PME grid + cut-off optimization with GPUs or PME nodes */
1881 /* Count the total cycles over the last steps */
1882 cycles_pmes += cycles;
1884 /* We can only switch cut-off at NS steps */
1885 if (step % ir->nstlist == 0)
1887 /* PME grid + cut-off optimization with GPUs or PME nodes */
1890 if (DDMASTER(cr->dd))
1892 /* PME node load is too high, start tuning */
1893 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1895 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1897 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1899 bPMETuneTry = FALSE;
1902 if (bPMETuneRunning)
1904 /* init_step might not be a multiple of nstlist,
1905 * but the first cycle is always skipped anyhow.
1908 pme_load_balance(pme_loadbal, cr,
1909 (bVerbose && MASTER(cr)) ? stderr : NULL,
1911 ir, state, cycles_pmes,
1912 fr->ic, fr->nbv, &fr->pmedata,
1915 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1916 fr->ewaldcoeff = fr->ic->ewaldcoeff;
1917 fr->rlist = fr->ic->rlist;
1918 fr->rlistlong = fr->ic->rlistlong;
1919 fr->rcoulomb = fr->ic->rcoulomb;
1920 fr->rvdw = fr->ic->rvdw;
1926 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1927 gs.set[eglsRESETCOUNTERS] != 0)
1929 /* Reset all the counters related to performance over the run */
1930 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, runtime,
1931 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1932 wcycle_set_reset_counters(wcycle, -1);
1933 if (!(cr->duty & DUTY_PME))
1935 /* Tell our PME node to reset its counters */
1936 gmx_pme_send_resetcounters(cr, step);
1938 /* Correct max_hours for the elapsed time */
1939 max_hours -= run_time/(60.0*60.0);
1940 bResetCountersHalfMaxH = FALSE;
1941 gs.set[eglsRESETCOUNTERS] = 0;
1945 /* End of main MD loop */
1949 runtime_end(runtime);
1951 if (bRerunMD && MASTER(cr))
1956 if (!(cr->duty & DUTY_PME))
1958 /* Tell the PME only node to finish */
1959 gmx_pme_send_finish(cr);
1964 if (ir->nstcalcenergy > 0 && !bRerunMD)
1966 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
1967 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1975 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1977 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)));
1978 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1981 if (pme_loadbal != NULL)
1983 pme_loadbal_done(pme_loadbal, cr, fplog,
1984 fr->nbv != NULL && fr->nbv->bUseGPU);
1987 if (shellfc && fplog)
1989 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1990 (nconverged*100.0)/step_rel);
1991 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1995 if (repl_ex_nst > 0 && MASTER(cr))
1997 print_replica_exchange_statistics(fplog, repl_ex);
2000 runtime->nsteps_done = step_rel;