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
53 #include "md_support.h"
54 #include "md_logging.h"
68 #include "domdec_network.h"
74 #include "compute_io.h"
76 #include "checkpoint.h"
77 #include "mtop_util.h"
78 #include "sighandler.h"
81 #include "pme_loadbal.h"
84 #include "types/nlistheuristics.h"
85 #include "types/iteratedconstraints.h"
86 #include "nbnxn_cuda_data_mgmt.h"
88 #include "gromacs/utility/gmxmpi.h"
89 #include "gromacs/fileio/confio.h"
90 #include "gromacs/fileio/trajectory_writing.h"
91 #include "gromacs/fileio/trnio.h"
92 #include "gromacs/fileio/trxio.h"
93 #include "gromacs/fileio/xtcio.h"
94 #include "gromacs/timing/wallcycle.h"
95 #include "gromacs/timing/walltime_accounting.h"
101 static void reset_all_counters(FILE *fplog, t_commrec *cr,
102 gmx_large_int_t step,
103 gmx_large_int_t *step_rel, t_inputrec *ir,
104 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
105 gmx_walltime_accounting_t walltime_accounting,
106 nbnxn_cuda_ptr_t cu_nbv)
108 char sbuf[STEPSTRSIZE];
110 /* Reset all the counters related to performance over the run */
111 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
112 gmx_step_str(step, sbuf));
116 nbnxn_cuda_reset_timings(cu_nbv);
119 wallcycle_stop(wcycle, ewcRUN);
120 wallcycle_reset_all(wcycle);
121 if (DOMAINDECOMP(cr))
123 reset_dd_statistics_counters(cr->dd);
126 ir->init_step += *step_rel;
127 ir->nsteps -= *step_rel;
129 wallcycle_start(wcycle, ewcRUN);
130 walltime_accounting_start(walltime_accounting);
131 print_date_and_time(fplog, cr->nodeid, "Restarted time", walltime_accounting);
134 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
135 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
137 gmx_vsite_t *vsite, gmx_constr_t constr,
138 int stepout, t_inputrec *ir,
139 gmx_mtop_t *top_global,
141 t_state *state_global,
143 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
144 gmx_edsam_t ed, t_forcerec *fr,
145 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
146 real cpt_period, real max_hours,
147 const char gmx_unused *deviceOptions,
149 gmx_walltime_accounting_t walltime_accounting)
152 gmx_large_int_t step, step_rel;
154 double t, t0, lam0[efptNR];
155 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
156 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
157 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
158 bBornRadii, bStartingFromCpt;
159 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
160 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
161 bForceUpdate = FALSE, bCPT;
162 gmx_bool bMasterState;
163 int force_flags, cglo_flags;
164 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
169 t_state *bufstate = NULL;
170 matrix *scale_tot, pcoupl_mu, M, ebox;
173 gmx_repl_ex_t repl_ex = NULL;
176 t_mdebin *mdebin = NULL;
177 t_state *state = NULL;
178 rvec *f_global = NULL;
179 gmx_enerdata_t *enerd;
181 gmx_global_stat_t gstat;
182 gmx_update_t upd = NULL;
183 t_graph *graph = NULL;
185 gmx_rng_t mcrng = NULL;
186 gmx_groups_t *groups;
187 gmx_ekindata_t *ekind, *ekind_save;
188 gmx_shellfc_t shellfc;
189 int count, nconverged = 0;
192 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
194 gmx_bool bResetCountersHalfMaxH = FALSE;
195 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
196 gmx_bool bUpdateDoLR;
201 real veta_save, scalevir, tracevir;
207 real saved_conserved_quantity = 0;
212 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
213 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
214 gmx_iterate_t iterate;
215 gmx_large_int_t multisim_nsteps = -1; /* number of steps to do before first multisim
216 simulation stops. If equal to zero, don't
217 communicate any more between multisims.*/
218 /* PME load balancing data for GPU kernels */
219 pme_load_balancing_t pme_loadbal = NULL;
221 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
224 /* Temporary addition for FAHCORE checkpointing */
228 /* Check for special mdrun options */
229 bRerunMD = (Flags & MD_RERUN);
230 bAppend = (Flags & MD_APPENDFILES);
231 if (Flags & MD_RESETCOUNTERSHALFWAY)
235 /* Signal to reset the counters half the simulation steps. */
236 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
238 /* Signal to reset the counters halfway the simulation time. */
239 bResetCountersHalfMaxH = (max_hours > 0);
242 /* md-vv uses averaged full step velocities for T-control
243 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
244 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
246 if (bVV) /* to store the initial velocities while computing virial */
248 snew(cbuf, top_global->natoms);
250 /* all the iteratative cases - only if there are constraints */
251 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
252 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
253 false in this step. The correct value, true or false,
254 is set at each step, as it depends on the frequency of temperature
255 and pressure control.*/
256 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
260 /* Since we don't know if the frames read are related in any way,
261 * rebuild the neighborlist at every step.
264 ir->nstcalcenergy = 1;
268 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
270 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
271 bGStatEveryStep = (nstglobalcomm == 1);
273 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
276 "To reduce the energy communication with nstlist = -1\n"
277 "the neighbor list validity should not be checked at every step,\n"
278 "this means that exact integration is not guaranteed.\n"
279 "The neighbor list validity is checked after:\n"
280 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
281 "In most cases this will result in exact integration.\n"
282 "This reduces the energy communication by a factor of 2 to 3.\n"
283 "If you want less energy communication, set nstlist > 3.\n\n");
290 groups = &top_global->groups;
293 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
294 &(state_global->fep_state), lam0,
295 nrnb, top_global, &upd,
296 nfile, fnm, &outf, &mdebin,
297 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
299 clear_mat(total_vir);
301 /* Energy terms and groups */
303 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
305 if (DOMAINDECOMP(cr))
311 snew(f, top_global->natoms);
314 /* Kinetic energy data */
316 init_ekindata(fplog, top_global, &(ir->opts), ekind);
317 /* needed for iteration of constraints */
319 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
320 /* Copy the cos acceleration to the groups struct */
321 ekind->cosacc.cos_accel = ir->cos_accel;
323 gstat = global_stat_init(ir);
326 /* Check for polarizable models and flexible constraints */
327 shellfc = init_shell_flexcon(fplog,
328 top_global, n_flexible_constraints(constr),
329 (ir->bContinuation ||
330 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
331 NULL : state_global->x);
335 #ifdef GMX_THREAD_MPI
336 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
338 set_deform_reference_box(upd,
339 deform_init_init_step_tpx,
340 deform_init_box_tpx);
341 #ifdef GMX_THREAD_MPI
342 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
347 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
348 if ((io > 2000) && MASTER(cr))
351 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
356 if (DOMAINDECOMP(cr))
358 top = dd_init_local_top(top_global);
361 dd_init_local_state(cr->dd, state_global, state);
363 if (DDMASTER(cr->dd) && ir->nstfout)
365 snew(f_global, state_global->natoms);
372 /* Initialize the particle decomposition and split the topology */
373 top = split_system(fplog, top_global, ir, cr);
375 pd_cg_range(cr, &fr->cg0, &fr->hcg);
376 pd_at_range(cr, &a0, &a1);
380 top = gmx_mtop_generate_local_top(top_global, ir);
383 a1 = top_global->natoms;
386 forcerec_set_excl_load(fr, top, cr);
388 state = partdec_init_local_state(cr, state_global);
391 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
395 set_vsite_top(vsite, top, mdatoms, cr);
398 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
400 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
405 make_local_shells(cr, mdatoms, shellfc);
408 setup_bonded_threading(fr, &top->idef);
410 if (ir->pull && PAR(cr))
412 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
416 if (DOMAINDECOMP(cr))
418 /* Distribute the charge groups over the nodes from the master node */
419 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
420 state_global, top_global, ir,
421 state, &f, mdatoms, top, fr,
422 vsite, shellfc, constr,
423 nrnb, wcycle, FALSE);
427 update_mdatoms(mdatoms, state->lambda[efptMASS]);
429 if (opt2bSet("-cpi", nfile, fnm))
431 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
435 bStateFromCP = FALSE;
440 init_expanded_ensemble(bStateFromCP,ir,&mcrng,&state->dfhist);
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 walltime_accounting_start(walltime_accounting);
676 print_date_and_time(fplog, cr->nodeid, "Started mdrun", walltime_accounting);
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)
814 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
819 update_annealing_target_temp(&(ir->opts), t);
824 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
826 for (i = 0; i < state_global->natoms; i++)
828 copy_rvec(rerun_fr.x[i], state_global->x[i]);
832 for (i = 0; i < state_global->natoms; i++)
834 copy_rvec(rerun_fr.v[i], state_global->v[i]);
839 for (i = 0; i < state_global->natoms; i++)
841 clear_rvec(state_global->v[i]);
845 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
846 " Ekin, temperature and pressure are incorrect,\n"
847 " the virial will be incorrect when constraints are present.\n"
849 bRerunWarnNoV = FALSE;
853 copy_mat(rerun_fr.box, state_global->box);
854 copy_mat(state_global->box, state->box);
856 if (vsite && (Flags & MD_RERUN_VSITE))
858 if (DOMAINDECOMP(cr))
860 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
864 /* Following is necessary because the graph may get out of sync
865 * with the coordinates if we only have every N'th coordinate set
867 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
868 shift_self(graph, state->box, state->x);
870 construct_vsites(vsite, state->x, ir->delta_t, state->v,
871 top->idef.iparams, top->idef.il,
872 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
875 unshift_self(graph, state->box, state->x);
880 /* Stop Center of Mass motion */
881 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
885 /* for rerun MD always do Neighbour Searching */
886 bNS = (bFirstStep || ir->nstlist != 0);
891 /* Determine whether or not to do Neighbour Searching and LR */
892 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
894 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
895 (ir->nstlist == -1 && nlh.nabnsb > 0));
897 if (bNS && ir->nstlist == -1)
899 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
903 /* check whether we should stop because another simulation has
907 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
908 (multisim_nsteps != ir->nsteps) )
915 "Stopping simulation %d because another one has finished\n",
919 gs.sig[eglsCHKPT] = 1;
924 /* < 0 means stop at next step, > 0 means stop at next NS step */
925 if ( (gs.set[eglsSTOPCOND] < 0) ||
926 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
931 /* Determine whether or not to update the Born radii if doing GB */
932 bBornRadii = bFirstStep;
933 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
938 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
939 do_verbose = bVerbose &&
940 (step % stepout == 0 || bFirstStep || bLastStep);
942 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
950 bMasterState = FALSE;
951 /* Correct the new box if it is too skewed */
952 if (DYNAMIC_BOX(*ir))
954 if (correct_box(fplog, step, state->box, graph))
959 if (DOMAINDECOMP(cr) && bMasterState)
961 dd_collect_state(cr->dd, state, state_global);
965 if (DOMAINDECOMP(cr))
967 /* Repartition the domain decomposition */
968 wallcycle_start(wcycle, ewcDOMDEC);
969 dd_partition_system(fplog, step, cr,
970 bMasterState, nstglobalcomm,
971 state_global, top_global, ir,
972 state, &f, mdatoms, top, fr,
973 vsite, shellfc, constr,
975 do_verbose && !bPMETuneRunning);
976 wallcycle_stop(wcycle, ewcDOMDEC);
977 /* If using an iterative integrator, reallocate space to match the decomposition */
981 if (MASTER(cr) && do_log)
983 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
986 if (ir->efep != efepNO)
988 update_mdatoms(mdatoms, state->lambda[efptMASS]);
991 if ((bRerunMD && rerun_fr.bV) || bExchanged)
994 /* We need the kinetic energy at minus the half step for determining
995 * the full step kinetic energy and possibly for T-coupling.*/
996 /* This may not be quite working correctly yet . . . . */
997 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
998 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
999 constr, NULL, FALSE, state->box,
1000 top_global, &bSumEkinhOld,
1001 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1003 clear_mat(force_vir);
1005 /* We write a checkpoint at this MD step when:
1006 * either at an NS step when we signalled through gs,
1007 * or at the last step (but not when we do not want confout),
1008 * but never at the first step or with rerun.
1010 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1011 (bLastStep && (Flags & MD_CONFOUT))) &&
1012 step > ir->init_step && !bRerunMD);
1015 gs.set[eglsCHKPT] = 0;
1018 /* Determine the energy and pressure:
1019 * at nstcalcenergy steps and at energy output steps (set below).
1021 if (EI_VV(ir->eI) && (!bInitStep))
1023 /* for vv, the first half of the integration actually corresponds
1024 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1025 but the virial needs to be calculated on both the current step and the 'next' step. Future
1026 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1028 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1029 bCalcVir = bCalcEner ||
1030 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1034 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1035 bCalcVir = bCalcEner ||
1036 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1039 /* Do we need global communication ? */
1040 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1041 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1042 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1044 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1046 if (do_ene || do_log)
1053 /* these CGLO_ options remain the same throughout the iteration */
1054 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1055 (bGStat ? CGLO_GSTAT : 0)
1058 force_flags = (GMX_FORCE_STATECHANGED |
1059 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1060 GMX_FORCE_ALLFORCES |
1062 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1063 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1064 (bDoFEP ? GMX_FORCE_DHDL : 0)
1069 if (do_per_step(step, ir->nstcalclr))
1071 force_flags |= GMX_FORCE_DO_LR;
1077 /* Now is the time to relax the shells */
1078 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1079 ir, bNS, force_flags,
1082 state, f, force_vir, mdatoms,
1083 nrnb, wcycle, graph, groups,
1084 shellfc, fr, bBornRadii, t, mu_tot,
1096 /* The coordinates (x) are shifted (to get whole molecules)
1098 * This is parallellized as well, and does communication too.
1099 * Check comments in sim_util.c
1101 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1102 state->box, state->x, &state->hist,
1103 f, force_vir, mdatoms, enerd, fcd,
1104 state->lambda, graph,
1105 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1106 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1109 if (bVV && !bStartingFromCpt && !bRerunMD)
1110 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1112 if (ir->eI == eiVV && bInitStep)
1114 /* if using velocity verlet with full time step Ekin,
1115 * take the first half step only to compute the
1116 * virial for the first step. From there,
1117 * revert back to the initial coordinates
1118 * so that the input is actually the initial step.
1120 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1124 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1125 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1128 /* If we are using twin-range interactions where the long-range component
1129 * is only evaluated every nstcalclr>1 steps, we should do a special update
1130 * step to combine the long-range forces on these steps.
1131 * For nstcalclr=1 this is not done, since the forces would have been added
1132 * directly to the short-range forces already.
1134 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1136 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1137 f, bUpdateDoLR, fr->f_twin, fcd,
1138 ekind, M, upd, bInitStep, etrtVELOCITY1,
1139 cr, nrnb, constr, &top->idef);
1141 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1143 gmx_iterate_init(&iterate, TRUE);
1145 /* for iterations, we save these vectors, as we will be self-consistently iterating
1148 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1150 /* save the state */
1151 if (iterate.bIterationActive)
1153 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1156 bFirstIterate = TRUE;
1157 while (bFirstIterate || iterate.bIterationActive)
1159 if (iterate.bIterationActive)
1161 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1162 if (bFirstIterate && bTrotter)
1164 /* The first time through, we need a decent first estimate
1165 of veta(t+dt) to compute the constraints. Do
1166 this by computing the box volume part of the
1167 trotter integration at this time. Nothing else
1168 should be changed by this routine here. If
1169 !(first time), we start with the previous value
1172 veta_save = state->veta;
1173 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1174 vetanew = state->veta;
1175 state->veta = veta_save;
1180 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1182 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1183 state, fr->bMolPBC, graph, f,
1184 &top->idef, shake_vir,
1185 cr, nrnb, wcycle, upd, constr,
1186 TRUE, bCalcVir, vetanew);
1190 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1196 /* Need to unshift here if a do_force has been
1197 called in the previous step */
1198 unshift_self(graph, state->box, state->x);
1201 /* if VV, compute the pressure and constraints */
1202 /* For VV2, we strictly only need this if using pressure
1203 * control, but we really would like to have accurate pressures
1205 * Think about ways around this in the future?
1206 * For now, keep this choice in comments.
1208 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1209 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1211 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1212 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1214 bSumEkinhOld = TRUE;
1216 /* for vv, the first half of the integration actually corresponds to the previous step.
1217 So we need information from the last step in the first half of the integration */
1218 if (bGStat || do_per_step(step-1, nstglobalcomm))
1220 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1221 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1222 constr, NULL, FALSE, state->box,
1223 top_global, &bSumEkinhOld,
1226 | (bTemp ? CGLO_TEMPERATURE : 0)
1227 | (bPres ? CGLO_PRESSURE : 0)
1228 | (bPres ? CGLO_CONSTRAINT : 0)
1229 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1230 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1233 /* explanation of above:
1234 a) We compute Ekin at the full time step
1235 if 1) we are using the AveVel Ekin, and it's not the
1236 initial step, or 2) if we are using AveEkin, but need the full
1237 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1238 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1239 EkinAveVel because it's needed for the pressure */
1241 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1246 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1247 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1254 /* We need the kinetic energy at minus the half step for determining
1255 * the full step kinetic energy and possibly for T-coupling.*/
1256 /* This may not be quite working correctly yet . . . . */
1257 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1258 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1259 constr, NULL, FALSE, state->box,
1260 top_global, &bSumEkinhOld,
1261 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1266 if (iterate.bIterationActive &&
1267 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1268 state->veta, &vetanew))
1272 bFirstIterate = FALSE;
1275 if (bTrotter && !bInitStep)
1277 copy_mat(shake_vir, state->svir_prev);
1278 copy_mat(force_vir, state->fvir_prev);
1279 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1281 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1282 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1283 enerd->term[F_EKIN] = trace(ekind->ekin);
1286 /* if it's the initial step, we performed this first step just to get the constraint virial */
1287 if (bInitStep && ir->eI == eiVV)
1289 copy_rvecn(cbuf, state->v, 0, state->natoms);
1293 /* MRS -- now done iterating -- compute the conserved quantity */
1296 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1299 last_ekin = enerd->term[F_EKIN];
1301 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1303 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1305 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1308 sum_dhdl(enerd, state->lambda, ir->fepvals);
1312 /* ######## END FIRST UPDATE STEP ############## */
1313 /* ######## If doing VV, we now have v(dt) ###### */
1316 /* perform extended ensemble sampling in lambda - we don't
1317 actually move to the new state before outputting
1318 statistics, but if performing simulated tempering, we
1319 do update the velocities and the tau_t. */
1321 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, mcrng, state->v, mdatoms);
1322 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1323 copy_df_history(&state_global->dfhist,&state->dfhist);
1326 /* Now we have the energies and forces corresponding to the
1327 * coordinates at time t. We must output all of this before
1330 do_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1331 ir, state, state_global, top_global, fr, upd,
1332 outf, mdebin, ekind, f, f_global,
1333 wcycle, mcrng, &nchkpt,
1334 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1337 /* kludge -- virial is lost with restart for NPT control. Must restart */
1338 if (bStartingFromCpt && bVV)
1340 copy_mat(state->svir_prev, shake_vir);
1341 copy_mat(state->fvir_prev, force_vir);
1344 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1346 /* Check whether everything is still allright */
1347 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1348 #ifdef GMX_THREAD_MPI
1353 /* this is just make gs.sig compatible with the hack
1354 of sending signals around by MPI_Reduce with together with
1356 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1358 gs.sig[eglsSTOPCOND] = 1;
1360 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1362 gs.sig[eglsSTOPCOND] = -1;
1364 /* < 0 means stop at next step, > 0 means stop at next NS step */
1368 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1369 gmx_get_signal_name(),
1370 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1374 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1375 gmx_get_signal_name(),
1376 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1378 handled_stop_condition = (int)gmx_get_stop_condition();
1380 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1381 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1382 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1384 /* Signal to terminate the run */
1385 gs.sig[eglsSTOPCOND] = 1;
1388 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1390 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1393 if (bResetCountersHalfMaxH && MASTER(cr) &&
1394 elapsed_time > max_hours*60.0*60.0*0.495)
1396 gs.sig[eglsRESETCOUNTERS] = 1;
1399 if (ir->nstlist == -1 && !bRerunMD)
1401 /* When bGStatEveryStep=FALSE, global_stat is only called
1402 * when we check the atom displacements, not at NS steps.
1403 * This means that also the bonded interaction count check is not
1404 * performed immediately after NS. Therefore a few MD steps could
1405 * be performed with missing interactions.
1406 * But wrong energies are never written to file,
1407 * since energies are only written after global_stat
1410 if (step >= nlh.step_nscheck)
1412 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1413 nlh.scale_tot, state->x);
1417 /* This is not necessarily true,
1418 * but step_nscheck is determined quite conservatively.
1424 /* In parallel we only have to check for checkpointing in steps
1425 * where we do global communication,
1426 * otherwise the other nodes don't know.
1428 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1431 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1432 gs.set[eglsCHKPT] == 0)
1434 gs.sig[eglsCHKPT] = 1;
1437 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1442 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1444 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1446 gmx_bool bIfRandomize;
1447 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1448 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1449 if (constr && bIfRandomize)
1451 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1452 state, fr->bMolPBC, graph, f,
1453 &top->idef, tmp_vir,
1454 cr, nrnb, wcycle, upd, constr,
1455 TRUE, bCalcVir, vetanew);
1460 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1462 gmx_iterate_init(&iterate, TRUE);
1463 /* for iterations, we save these vectors, as we will be redoing the calculations */
1464 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1467 bFirstIterate = TRUE;
1468 while (bFirstIterate || iterate.bIterationActive)
1470 /* We now restore these vectors to redo the calculation with improved extended variables */
1471 if (iterate.bIterationActive)
1473 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1476 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1477 so scroll down for that logic */
1479 /* ######### START SECOND UPDATE STEP ################# */
1480 /* Box is changed in update() when we do pressure coupling,
1481 * but we should still use the old box for energy corrections and when
1482 * writing it to the energy file, so it matches the trajectory files for
1483 * the same timestep above. Make a copy in a separate array.
1485 copy_mat(state->box, lastbox);
1490 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1492 wallcycle_start(wcycle, ewcUPDATE);
1493 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1496 if (iterate.bIterationActive)
1504 /* we use a new value of scalevir to converge the iterations faster */
1505 scalevir = tracevir/trace(shake_vir);
1507 msmul(shake_vir, scalevir, shake_vir);
1508 m_add(force_vir, shake_vir, total_vir);
1509 clear_mat(shake_vir);
1511 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1512 /* We can only do Berendsen coupling after we have summed
1513 * the kinetic energy or virial. Since the happens
1514 * in global_state after update, we should only do it at
1515 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1520 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1521 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1526 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1528 /* velocity half-step update */
1529 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1530 bUpdateDoLR, fr->f_twin, fcd,
1531 ekind, M, upd, FALSE, etrtVELOCITY2,
1532 cr, nrnb, constr, &top->idef);
1535 /* Above, initialize just copies ekinh into ekin,
1536 * it doesn't copy position (for VV),
1537 * and entire integrator for MD.
1540 if (ir->eI == eiVVAK)
1542 copy_rvecn(state->x, cbuf, 0, state->natoms);
1544 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1546 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1547 bUpdateDoLR, fr->f_twin, fcd,
1548 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1549 wallcycle_stop(wcycle, ewcUPDATE);
1551 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1552 fr->bMolPBC, graph, f,
1553 &top->idef, shake_vir,
1554 cr, nrnb, wcycle, upd, constr,
1555 FALSE, bCalcVir, state->veta);
1557 if (ir->eI == eiVVAK)
1559 /* erase F_EKIN and F_TEMP here? */
1560 /* just compute the kinetic energy at the half step to perform a trotter step */
1561 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1562 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1563 constr, NULL, FALSE, lastbox,
1564 top_global, &bSumEkinhOld,
1565 cglo_flags | CGLO_TEMPERATURE
1567 wallcycle_start(wcycle, ewcUPDATE);
1568 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1569 /* now we know the scaling, we can compute the positions again again */
1570 copy_rvecn(cbuf, state->x, 0, state->natoms);
1572 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1574 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1575 bUpdateDoLR, fr->f_twin, fcd,
1576 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1577 wallcycle_stop(wcycle, ewcUPDATE);
1579 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1580 /* are the small terms in the shake_vir here due
1581 * to numerical errors, or are they important
1582 * physically? I'm thinking they are just errors, but not completely sure.
1583 * For now, will call without actually constraining, constr=NULL*/
1584 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1585 state, fr->bMolPBC, graph, f,
1586 &top->idef, tmp_vir,
1587 cr, nrnb, wcycle, upd, NULL,
1593 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1596 if (fr->bSepDVDL && fplog && do_log)
1598 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1602 /* this factor or 2 correction is necessary
1603 because half of the constraint force is removed
1604 in the vv step, so we have to double it. See
1605 the Redmine issue #1255. It is not yet clear
1606 if the factor of 2 is exact, or just a very
1607 good approximation, and this will be
1608 investigated. The next step is to see if this
1609 can be done adding a dhdl contribution from the
1610 rattle step, but this is somewhat more
1611 complicated with the current code. Will be
1612 investigated, hopefully for 4.6.3. However,
1613 this current solution is much better than
1614 having it completely wrong.
1616 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1620 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1625 /* Need to unshift here */
1626 unshift_self(graph, state->box, state->x);
1631 wallcycle_start(wcycle, ewcVSITECONSTR);
1634 shift_self(graph, state->box, state->x);
1636 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1637 top->idef.iparams, top->idef.il,
1638 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1642 unshift_self(graph, state->box, state->x);
1644 wallcycle_stop(wcycle, ewcVSITECONSTR);
1647 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1648 /* With Leap-Frog we can skip compute_globals at
1649 * non-communication steps, but we need to calculate
1650 * the kinetic energy one step before communication.
1652 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1654 if (ir->nstlist == -1 && bFirstIterate)
1656 gs.sig[eglsNABNSB] = nlh.nabnsb;
1658 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1659 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1661 bFirstIterate ? &gs : NULL,
1662 (step_rel % gs.nstms == 0) &&
1663 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1665 top_global, &bSumEkinhOld,
1667 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1668 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1669 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1670 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1671 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1672 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1675 if (ir->nstlist == -1 && bFirstIterate)
1677 nlh.nabnsb = gs.set[eglsNABNSB];
1678 gs.set[eglsNABNSB] = 0;
1681 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1682 /* ############# END CALC EKIN AND PRESSURE ################# */
1684 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1685 the virial that should probably be addressed eventually. state->veta has better properies,
1686 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1687 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1689 if (iterate.bIterationActive &&
1690 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1691 trace(shake_vir), &tracevir))
1695 bFirstIterate = FALSE;
1698 if (!bVV || bRerunMD)
1700 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1701 sum_dhdl(enerd, state->lambda, ir->fepvals);
1703 update_box(fplog, step, ir, mdatoms, state, f,
1704 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1706 /* ################# END UPDATE STEP 2 ################# */
1707 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1709 /* The coordinates (x) were unshifted in update */
1712 /* We will not sum ekinh_old,
1713 * so signal that we still have to do it.
1715 bSumEkinhOld = TRUE;
1718 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1720 /* use the directly determined last velocity, not actually the averaged half steps */
1721 if (bTrotter && ir->eI == eiVV)
1723 enerd->term[F_EKIN] = last_ekin;
1725 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1729 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1733 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1735 /* ######### END PREPARING EDR OUTPUT ########### */
1740 gmx_bool do_dr, do_or;
1742 if (fplog && do_log && bDoExpanded)
1744 /* only needed if doing expanded ensemble */
1745 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1746 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1748 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1752 upd_mdebin(mdebin, bDoDHDL, TRUE,
1753 t, mdatoms->tmass, enerd, state,
1754 ir->fepvals, ir->expandedvals, lastbox,
1755 shake_vir, force_vir, total_vir, pres,
1756 ekind, mu_tot, constr);
1760 upd_mdebin_step(mdebin);
1763 do_dr = do_per_step(step, ir->nstdisreout);
1764 do_or = do_per_step(step, ir->nstorireout);
1766 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
1768 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1770 if (ir->ePull != epullNO)
1772 pull_print_output(ir->pull, step, t);
1775 if (do_per_step(step, ir->nstlog))
1777 if (fflush(fplog) != 0)
1779 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1785 /* Have to do this part _after_ outputting the logfile and the edr file */
1786 /* Gets written into the state at the beginning of next loop*/
1787 state->fep_state = lamnew;
1789 /* Print the remaining wall clock time for the run */
1790 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1794 fprintf(stderr, "\n");
1796 print_time(stderr, walltime_accounting, step, ir, cr);
1799 /* Replica exchange */
1801 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1802 do_per_step(step, repl_ex_nst))
1804 bExchanged = replica_exchange(fplog, cr, repl_ex,
1805 state_global, enerd,
1808 if (bExchanged && DOMAINDECOMP(cr))
1810 dd_partition_system(fplog, step, cr, TRUE, 1,
1811 state_global, top_global, ir,
1812 state, &f, mdatoms, top, fr,
1813 vsite, shellfc, constr,
1814 nrnb, wcycle, FALSE);
1820 bStartingFromCpt = FALSE;
1822 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1823 /* With all integrators, except VV, we need to retain the pressure
1824 * at the current step for coupling at the next step.
1826 if ((state->flags & (1<<estPRES_PREV)) &&
1828 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1830 /* Store the pressure in t_state for pressure coupling
1831 * at the next MD step.
1833 copy_mat(pres, state->pres_prev);
1836 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1838 if ( (membed != NULL) && (!bLastStep) )
1840 rescale_membed(step_rel, membed, state_global->x);
1847 /* read next frame from input trajectory */
1848 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1853 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1857 if (!bRerunMD || !rerun_fr.bStep)
1859 /* increase the MD step number */
1864 cycles = wallcycle_stop(wcycle, ewcSTEP);
1865 if (DOMAINDECOMP(cr) && wcycle)
1867 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1870 if (bPMETuneRunning || bPMETuneTry)
1872 /* PME grid + cut-off optimization with GPUs or PME nodes */
1874 /* Count the total cycles over the last steps */
1875 cycles_pmes += cycles;
1877 /* We can only switch cut-off at NS steps */
1878 if (step % ir->nstlist == 0)
1880 /* PME grid + cut-off optimization with GPUs or PME nodes */
1883 if (DDMASTER(cr->dd))
1885 /* PME node load is too high, start tuning */
1886 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1888 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1890 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1892 bPMETuneTry = FALSE;
1895 if (bPMETuneRunning)
1897 /* init_step might not be a multiple of nstlist,
1898 * but the first cycle is always skipped anyhow.
1901 pme_load_balance(pme_loadbal, cr,
1902 (bVerbose && MASTER(cr)) ? stderr : NULL,
1904 ir, state, cycles_pmes,
1905 fr->ic, fr->nbv, &fr->pmedata,
1908 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1909 fr->ewaldcoeff = fr->ic->ewaldcoeff;
1910 fr->rlist = fr->ic->rlist;
1911 fr->rlistlong = fr->ic->rlistlong;
1912 fr->rcoulomb = fr->ic->rcoulomb;
1913 fr->rvdw = fr->ic->rvdw;
1919 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1920 gs.set[eglsRESETCOUNTERS] != 0)
1922 /* Reset all the counters related to performance over the run */
1923 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1924 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1925 wcycle_set_reset_counters(wcycle, -1);
1926 if (!(cr->duty & DUTY_PME))
1928 /* Tell our PME node to reset its counters */
1929 gmx_pme_send_resetcounters(cr, step);
1931 /* Correct max_hours for the elapsed time */
1932 max_hours -= elapsed_time/(60.0*60.0);
1933 bResetCountersHalfMaxH = FALSE;
1934 gs.set[eglsRESETCOUNTERS] = 0;
1938 /* End of main MD loop */
1941 /* Stop measuring walltime */
1942 walltime_accounting_end(walltime_accounting);
1944 if (bRerunMD && MASTER(cr))
1949 if (!(cr->duty & DUTY_PME))
1951 /* Tell the PME only node to finish */
1952 gmx_pme_send_finish(cr);
1957 if (ir->nstcalcenergy > 0 && !bRerunMD)
1959 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
1960 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1968 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1970 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)));
1971 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1974 if (pme_loadbal != NULL)
1976 pme_loadbal_done(pme_loadbal, cr, fplog,
1977 fr->nbv != NULL && fr->nbv->bUseGPU);
1980 if (shellfc && fplog)
1982 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1983 (nconverged*100.0)/step_rel);
1984 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1988 if (repl_ex_nst > 0 && MASTER(cr))
1990 print_replica_exchange_statistics(fplog, repl_ex);
1993 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);