2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
54 #include "md_support.h"
55 #include "md_logging.h"
69 #include "domdec_network.h"
75 #include "compute_io.h"
77 #include "checkpoint.h"
78 #include "mtop_util.h"
79 #include "sighandler.h"
82 #include "pme_loadbal.h"
85 #include "types/nlistheuristics.h"
86 #include "types/iteratedconstraints.h"
87 #include "nbnxn_cuda_data_mgmt.h"
89 #include "gromacs/utility/gmxmpi.h"
90 #include "gromacs/fileio/confio.h"
91 #include "gromacs/fileio/trajectory_writing.h"
92 #include "gromacs/fileio/trnio.h"
93 #include "gromacs/fileio/trxio.h"
94 #include "gromacs/fileio/xtcio.h"
95 #include "gromacs/timing/wallcycle.h"
96 #include "gromacs/timing/walltime_accounting.h"
102 static void reset_all_counters(FILE *fplog, t_commrec *cr,
104 gmx_int64_t *step_rel, t_inputrec *ir,
105 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
106 gmx_walltime_accounting_t walltime_accounting,
107 nbnxn_cuda_ptr_t cu_nbv)
109 char sbuf[STEPSTRSIZE];
111 /* Reset all the counters related to performance over the run */
112 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
113 gmx_step_str(step, sbuf));
117 nbnxn_cuda_reset_timings(cu_nbv);
120 wallcycle_stop(wcycle, ewcRUN);
121 wallcycle_reset_all(wcycle);
122 if (DOMAINDECOMP(cr))
124 reset_dd_statistics_counters(cr->dd);
127 ir->init_step += *step_rel;
128 ir->nsteps -= *step_rel;
130 wallcycle_start(wcycle, ewcRUN);
131 walltime_accounting_start(walltime_accounting);
132 print_date_and_time(fplog, cr->nodeid, "Restarted time", walltime_accounting);
135 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
136 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
138 gmx_vsite_t *vsite, gmx_constr_t constr,
139 int stepout, t_inputrec *ir,
140 gmx_mtop_t *top_global,
142 t_state *state_global,
144 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
145 gmx_edsam_t ed, t_forcerec *fr,
146 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
147 real cpt_period, real max_hours,
148 const char gmx_unused *deviceOptions,
150 gmx_walltime_accounting_t walltime_accounting)
153 gmx_int64_t step, step_rel;
155 double t, t0, lam0[efptNR];
156 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
157 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
158 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
159 bBornRadii, bStartingFromCpt;
160 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
161 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
162 bForceUpdate = FALSE, bCPT;
163 gmx_bool bMasterState;
164 int force_flags, cglo_flags;
165 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
170 t_state *bufstate = NULL;
171 matrix *scale_tot, pcoupl_mu, M, ebox;
174 gmx_repl_ex_t repl_ex = NULL;
177 t_mdebin *mdebin = NULL;
178 t_state *state = NULL;
179 rvec *f_global = NULL;
180 gmx_enerdata_t *enerd;
182 gmx_global_stat_t gstat;
183 gmx_update_t upd = NULL;
184 t_graph *graph = NULL;
186 gmx_rng_t mcrng = NULL;
187 gmx_groups_t *groups;
188 gmx_ekindata_t *ekind, *ekind_save;
189 gmx_shellfc_t shellfc;
190 int count, nconverged = 0;
193 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
195 gmx_bool bResetCountersHalfMaxH = FALSE;
196 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
197 gmx_bool bUpdateDoLR;
202 real veta_save, scalevir, tracevir;
208 real saved_conserved_quantity = 0;
213 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
214 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
215 gmx_iterate_t iterate;
216 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
217 simulation stops. If equal to zero, don't
218 communicate any more between multisims.*/
219 /* PME load balancing data for GPU kernels */
220 pme_load_balancing_t pme_loadbal = NULL;
222 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
225 /* Temporary addition for FAHCORE checkpointing */
229 /* Check for special mdrun options */
230 bRerunMD = (Flags & MD_RERUN);
231 bAppend = (Flags & MD_APPENDFILES);
232 if (Flags & MD_RESETCOUNTERSHALFWAY)
236 /* Signal to reset the counters half the simulation steps. */
237 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
239 /* Signal to reset the counters halfway the simulation time. */
240 bResetCountersHalfMaxH = (max_hours > 0);
243 /* md-vv uses averaged full step velocities for T-control
244 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
245 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
247 if (bVV) /* to store the initial velocities while computing virial */
249 snew(cbuf, top_global->natoms);
251 /* all the iteratative cases - only if there are constraints */
252 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
253 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
254 false in this step. The correct value, true or false,
255 is set at each step, as it depends on the frequency of temperature
256 and pressure control.*/
257 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
261 /* Since we don't know if the frames read are related in any way,
262 * rebuild the neighborlist at every step.
265 ir->nstcalcenergy = 1;
269 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
271 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
272 bGStatEveryStep = (nstglobalcomm == 1);
274 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
277 "To reduce the energy communication with nstlist = -1\n"
278 "the neighbor list validity should not be checked at every step,\n"
279 "this means that exact integration is not guaranteed.\n"
280 "The neighbor list validity is checked after:\n"
281 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
282 "In most cases this will result in exact integration.\n"
283 "This reduces the energy communication by a factor of 2 to 3.\n"
284 "If you want less energy communication, set nstlist > 3.\n\n");
291 groups = &top_global->groups;
294 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
295 &(state_global->fep_state), lam0,
296 nrnb, top_global, &upd,
297 nfile, fnm, &outf, &mdebin,
298 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
300 clear_mat(total_vir);
302 /* Energy terms and groups */
304 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
306 if (DOMAINDECOMP(cr))
312 snew(f, top_global->natoms);
315 /* Kinetic energy data */
317 init_ekindata(fplog, top_global, &(ir->opts), ekind);
318 /* needed for iteration of constraints */
320 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
321 /* Copy the cos acceleration to the groups struct */
322 ekind->cosacc.cos_accel = ir->cos_accel;
324 gstat = global_stat_init(ir);
327 /* Check for polarizable models and flexible constraints */
328 shellfc = init_shell_flexcon(fplog,
329 top_global, n_flexible_constraints(constr),
330 (ir->bContinuation ||
331 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
332 NULL : state_global->x);
336 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
337 set_deform_reference_box(upd,
338 deform_init_init_step_tpx,
339 deform_init_box_tpx);
340 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
344 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
345 if ((io > 2000) && MASTER(cr))
348 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
353 if (DOMAINDECOMP(cr))
355 top = dd_init_local_top(top_global);
358 dd_init_local_state(cr->dd, state_global, state);
360 if (DDMASTER(cr->dd) && ir->nstfout)
362 snew(f_global, state_global->natoms);
369 /* Initialize the particle decomposition and split the topology */
370 top = split_system(fplog, top_global, ir, cr);
372 pd_cg_range(cr, &fr->cg0, &fr->hcg);
373 pd_at_range(cr, &a0, &a1);
377 top = gmx_mtop_generate_local_top(top_global, ir);
380 a1 = top_global->natoms;
383 forcerec_set_excl_load(fr, top, cr);
385 state = partdec_init_local_state(cr, state_global);
388 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
392 set_vsite_top(vsite, top, mdatoms, cr);
395 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
397 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
402 make_local_shells(cr, mdatoms, shellfc);
405 setup_bonded_threading(fr, &top->idef);
407 if (ir->pull && PAR(cr))
409 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
413 if (DOMAINDECOMP(cr))
415 /* Distribute the charge groups over the nodes from the master node */
416 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
417 state_global, top_global, ir,
418 state, &f, mdatoms, top, fr,
419 vsite, shellfc, constr,
420 nrnb, wcycle, FALSE);
424 update_mdatoms(mdatoms, state->lambda[efptMASS]);
426 if (opt2bSet("-cpi", nfile, fnm))
428 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
432 bStateFromCP = FALSE;
437 init_expanded_ensemble(bStateFromCP, ir, &mcrng, &state->dfhist);
444 /* Update mdebin with energy history if appending to output files */
445 if (Flags & MD_APPENDFILES)
447 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
451 /* We might have read an energy history from checkpoint,
452 * free the allocated memory and reset the counts.
454 done_energyhistory(&state_global->enerhist);
455 init_energyhistory(&state_global->enerhist);
458 /* Set the initial energy history in state by updating once */
459 update_energyhistory(&state_global->enerhist, mdebin);
462 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
464 /* Set the random state if we read a checkpoint file */
465 set_stochd_state(upd, state);
468 if (state->flags & (1<<estMC_RNG))
470 set_mc_state(mcrng, state);
473 /* Initialize constraints */
476 if (!DOMAINDECOMP(cr))
478 set_constraints(constr, top, ir, mdatoms, cr);
484 /* We need to be sure replica exchange can only occur
485 * when the energies are current */
486 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
487 "repl_ex_nst", &repl_ex_nst);
488 /* This check needs to happen before inter-simulation
489 * signals are initialized, too */
491 if (repl_ex_nst > 0 && MASTER(cr))
493 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
494 repl_ex_nst, repl_ex_nex, repl_ex_seed);
497 /* PME tuning is only supported with GPUs or PME nodes and not with rerun or LJ-PME.
498 * With perturbed charges with soft-core we should not change the cut-off.
500 if ((Flags & MD_TUNEPME) &&
501 EEL_PME(fr->eeltype) &&
502 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
503 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
504 !bRerunMD && !EVDW_PME(fr->vdwtype))
506 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
508 if (cr->duty & DUTY_PME)
510 /* Start tuning right away, as we can't measure the load */
511 bPMETuneRunning = TRUE;
515 /* Separate PME nodes, we can measure the PP/PME load balance */
520 if (!ir->bContinuation && !bRerunMD)
522 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
524 /* Set the velocities of frozen particles to zero */
525 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
527 for (m = 0; m < DIM; m++)
529 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
539 /* Constrain the initial coordinates and velocities */
540 do_constrain_first(fplog, constr, ir, mdatoms, state,
545 /* Construct the virtual sites for the initial configuration */
546 construct_vsites(vsite, state->x, ir->delta_t, NULL,
547 top->idef.iparams, top->idef.il,
548 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
554 /* set free energy calculation frequency as the minimum
555 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
556 nstfep = ir->fepvals->nstdhdl;
559 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
563 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
566 /* I'm assuming we need global communication the first time! MRS */
567 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
568 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
569 | (bVV ? CGLO_PRESSURE : 0)
570 | (bVV ? CGLO_CONSTRAINT : 0)
571 | (bRerunMD ? CGLO_RERUNMD : 0)
572 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
574 bSumEkinhOld = FALSE;
575 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
576 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
577 constr, NULL, FALSE, state->box,
578 top_global, &bSumEkinhOld, cglo_flags);
579 if (ir->eI == eiVVAK)
581 /* a second call to get the half step temperature initialized as well */
582 /* we do the same call as above, but turn the pressure off -- internally to
583 compute_globals, this is recognized as a velocity verlet half-step
584 kinetic energy calculation. This minimized excess variables, but
585 perhaps loses some logic?*/
587 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
588 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
589 constr, NULL, FALSE, state->box,
590 top_global, &bSumEkinhOld,
591 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
594 /* Calculate the initial half step temperature, and save the ekinh_old */
595 if (!(Flags & MD_STARTFROMCPT))
597 for (i = 0; (i < ir->opts.ngtc); i++)
599 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
604 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
605 and there is no previous step */
608 /* if using an iterative algorithm, we need to create a working directory for the state. */
611 bufstate = init_bufstate(state);
614 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
615 temperature control */
616 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
620 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
623 "RMS relative constraint deviation after constraining: %.2e\n",
624 constr_rmsd(constr, FALSE));
626 if (EI_STATE_VELOCITY(ir->eI))
628 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
632 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
633 " input trajectory '%s'\n\n",
634 *(top_global->name), opt2fn("-rerun", nfile, fnm));
637 fprintf(stderr, "Calculated time to finish depends on nsteps from "
638 "run input file,\nwhich may not correspond to the time "
639 "needed to process input trajectory.\n\n");
645 fprintf(stderr, "starting mdrun '%s'\n",
646 *(top_global->name));
649 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
653 sprintf(tbuf, "%s", "infinite");
655 if (ir->init_step > 0)
657 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
658 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
659 gmx_step_str(ir->init_step, sbuf2),
660 ir->init_step*ir->delta_t);
664 fprintf(stderr, "%s steps, %s ps.\n",
665 gmx_step_str(ir->nsteps, sbuf), tbuf);
668 fprintf(fplog, "\n");
671 /* Set and write start time */
672 walltime_accounting_start(walltime_accounting);
673 print_date_and_time(fplog, cr->nodeid, "Started mdrun", walltime_accounting);
674 wallcycle_start(wcycle, ewcRUN);
677 fprintf(fplog, "\n");
680 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
682 chkpt_ret = fcCheckPointParallel( cr->nodeid,
686 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
691 /***********************************************************
695 ************************************************************/
697 /* if rerunMD then read coordinates and velocities from input trajectory */
700 if (getenv("GMX_FORCE_UPDATE"))
708 bNotLastFrame = read_first_frame(oenv, &status,
709 opt2fn("-rerun", nfile, fnm),
710 &rerun_fr, TRX_NEED_X | TRX_READ_V);
711 if (rerun_fr.natoms != top_global->natoms)
714 "Number of atoms in trajectory (%d) does not match the "
715 "run input file (%d)\n",
716 rerun_fr.natoms, top_global->natoms);
718 if (ir->ePBC != epbcNONE)
722 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);
724 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
726 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
733 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
736 if (ir->ePBC != epbcNONE)
738 /* Set the shift vectors.
739 * Necessary here when have a static box different from the tpr box.
741 calc_shifts(rerun_fr.box, fr->shift_vec);
745 /* loop over MD steps or if rerunMD to end of input trajectory */
747 /* Skip the first Nose-Hoover integration when we get the state from tpx */
748 bStateFromTPX = !bStateFromCP;
749 bInitStep = bFirstStep && (bStateFromTPX || bVV);
750 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
752 bSumEkinhOld = FALSE;
755 init_global_signals(&gs, cr, ir, repl_ex_nst);
757 step = ir->init_step;
760 if (ir->nstlist == -1)
762 init_nlistheuristics(&nlh, bGStatEveryStep, step);
765 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
767 /* check how many steps are left in other sims */
768 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
772 /* and stop now if we should */
773 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
774 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
775 while (!bLastStep || (bRerunMD && bNotLastFrame))
778 wallcycle_start(wcycle, ewcSTEP);
784 step = rerun_fr.step;
785 step_rel = step - ir->init_step;
798 bLastStep = (step_rel == ir->nsteps);
799 t = t0 + step*ir->delta_t;
802 if (ir->efep != efepNO || ir->bSimTemp)
804 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
805 requiring different logic. */
807 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
808 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
809 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
810 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
811 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
816 update_annealing_target_temp(&(ir->opts), t);
821 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
823 for (i = 0; i < state_global->natoms; i++)
825 copy_rvec(rerun_fr.x[i], state_global->x[i]);
829 for (i = 0; i < state_global->natoms; i++)
831 copy_rvec(rerun_fr.v[i], state_global->v[i]);
836 for (i = 0; i < state_global->natoms; i++)
838 clear_rvec(state_global->v[i]);
842 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
843 " Ekin, temperature and pressure are incorrect,\n"
844 " the virial will be incorrect when constraints are present.\n"
846 bRerunWarnNoV = FALSE;
850 copy_mat(rerun_fr.box, state_global->box);
851 copy_mat(state_global->box, state->box);
853 if (vsite && (Flags & MD_RERUN_VSITE))
855 if (DOMAINDECOMP(cr))
857 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
861 /* Following is necessary because the graph may get out of sync
862 * with the coordinates if we only have every N'th coordinate set
864 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
865 shift_self(graph, state->box, state->x);
867 construct_vsites(vsite, state->x, ir->delta_t, state->v,
868 top->idef.iparams, top->idef.il,
869 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
872 unshift_self(graph, state->box, state->x);
877 /* Stop Center of Mass motion */
878 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
882 /* for rerun MD always do Neighbour Searching */
883 bNS = (bFirstStep || ir->nstlist != 0);
888 /* Determine whether or not to do Neighbour Searching and LR */
889 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
891 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
892 (ir->nstlist == -1 && nlh.nabnsb > 0));
894 if (bNS && ir->nstlist == -1)
896 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
900 /* check whether we should stop because another simulation has
904 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
905 (multisim_nsteps != ir->nsteps) )
912 "Stopping simulation %d because another one has finished\n",
916 gs.sig[eglsCHKPT] = 1;
921 /* < 0 means stop at next step, > 0 means stop at next NS step */
922 if ( (gs.set[eglsSTOPCOND] < 0) ||
923 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
928 /* Determine whether or not to update the Born radii if doing GB */
929 bBornRadii = bFirstStep;
930 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
935 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
936 do_verbose = bVerbose &&
937 (step % stepout == 0 || bFirstStep || bLastStep);
939 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
947 bMasterState = FALSE;
948 /* Correct the new box if it is too skewed */
949 if (DYNAMIC_BOX(*ir))
951 if (correct_box(fplog, step, state->box, graph))
956 if (DOMAINDECOMP(cr) && bMasterState)
958 dd_collect_state(cr->dd, state, state_global);
962 if (DOMAINDECOMP(cr))
964 /* Repartition the domain decomposition */
965 wallcycle_start(wcycle, ewcDOMDEC);
966 dd_partition_system(fplog, step, cr,
967 bMasterState, nstglobalcomm,
968 state_global, top_global, ir,
969 state, &f, mdatoms, top, fr,
970 vsite, shellfc, constr,
972 do_verbose && !bPMETuneRunning);
973 wallcycle_stop(wcycle, ewcDOMDEC);
974 /* If using an iterative integrator, reallocate space to match the decomposition */
978 if (MASTER(cr) && do_log)
980 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
983 if (ir->efep != efepNO)
985 update_mdatoms(mdatoms, state->lambda[efptMASS]);
988 if ((bRerunMD && rerun_fr.bV) || bExchanged)
991 /* We need the kinetic energy at minus the half step for determining
992 * the full step kinetic energy and possibly for T-coupling.*/
993 /* This may not be quite working correctly yet . . . . */
994 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
995 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
996 constr, NULL, FALSE, state->box,
997 top_global, &bSumEkinhOld,
998 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1000 clear_mat(force_vir);
1002 /* We write a checkpoint at this MD step when:
1003 * either at an NS step when we signalled through gs,
1004 * or at the last step (but not when we do not want confout),
1005 * but never at the first step or with rerun.
1007 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1008 (bLastStep && (Flags & MD_CONFOUT))) &&
1009 step > ir->init_step && !bRerunMD);
1012 gs.set[eglsCHKPT] = 0;
1015 /* Determine the energy and pressure:
1016 * at nstcalcenergy steps and at energy output steps (set below).
1018 if (EI_VV(ir->eI) && (!bInitStep))
1020 /* for vv, the first half of the integration actually corresponds
1021 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1022 but the virial needs to be calculated on both the current step and the 'next' step. Future
1023 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1025 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1026 bCalcVir = bCalcEner ||
1027 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1031 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1032 bCalcVir = bCalcEner ||
1033 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1036 /* Do we need global communication ? */
1037 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1038 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1039 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1041 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1043 if (do_ene || do_log)
1050 /* these CGLO_ options remain the same throughout the iteration */
1051 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1052 (bGStat ? CGLO_GSTAT : 0)
1055 force_flags = (GMX_FORCE_STATECHANGED |
1056 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1057 GMX_FORCE_ALLFORCES |
1059 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1060 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1061 (bDoFEP ? GMX_FORCE_DHDL : 0)
1066 if (do_per_step(step, ir->nstcalclr))
1068 force_flags |= GMX_FORCE_DO_LR;
1074 /* Now is the time to relax the shells */
1075 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1076 ir, bNS, force_flags,
1079 state, f, force_vir, mdatoms,
1080 nrnb, wcycle, graph, groups,
1081 shellfc, fr, bBornRadii, t, mu_tot,
1093 /* The coordinates (x) are shifted (to get whole molecules)
1095 * This is parallellized as well, and does communication too.
1096 * Check comments in sim_util.c
1098 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1099 state->box, state->x, &state->hist,
1100 f, force_vir, mdatoms, enerd, fcd,
1101 state->lambda, graph,
1102 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1103 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1106 if (bVV && !bStartingFromCpt && !bRerunMD)
1107 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1109 if (ir->eI == eiVV && bInitStep)
1111 /* if using velocity verlet with full time step Ekin,
1112 * take the first half step only to compute the
1113 * virial for the first step. From there,
1114 * revert back to the initial coordinates
1115 * so that the input is actually the initial step.
1117 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1121 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1122 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1125 /* If we are using twin-range interactions where the long-range component
1126 * is only evaluated every nstcalclr>1 steps, we should do a special update
1127 * step to combine the long-range forces on these steps.
1128 * For nstcalclr=1 this is not done, since the forces would have been added
1129 * directly to the short-range forces already.
1131 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1133 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1134 f, bUpdateDoLR, fr->f_twin, fcd,
1135 ekind, M, upd, bInitStep, etrtVELOCITY1,
1136 cr, nrnb, constr, &top->idef);
1138 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1140 gmx_iterate_init(&iterate, TRUE);
1142 /* for iterations, we save these vectors, as we will be self-consistently iterating
1145 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1147 /* save the state */
1148 if (iterate.bIterationActive)
1150 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1153 bFirstIterate = TRUE;
1154 while (bFirstIterate || iterate.bIterationActive)
1156 if (iterate.bIterationActive)
1158 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1159 if (bFirstIterate && bTrotter)
1161 /* The first time through, we need a decent first estimate
1162 of veta(t+dt) to compute the constraints. Do
1163 this by computing the box volume part of the
1164 trotter integration at this time. Nothing else
1165 should be changed by this routine here. If
1166 !(first time), we start with the previous value
1169 veta_save = state->veta;
1170 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1171 vetanew = state->veta;
1172 state->veta = veta_save;
1177 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1179 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1180 state, fr->bMolPBC, graph, f,
1181 &top->idef, shake_vir,
1182 cr, nrnb, wcycle, upd, constr,
1183 TRUE, bCalcVir, vetanew);
1187 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1193 /* Need to unshift here if a do_force has been
1194 called in the previous step */
1195 unshift_self(graph, state->box, state->x);
1198 /* if VV, compute the pressure and constraints */
1199 /* For VV2, we strictly only need this if using pressure
1200 * control, but we really would like to have accurate pressures
1202 * Think about ways around this in the future?
1203 * For now, keep this choice in comments.
1205 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1206 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1208 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1209 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1211 bSumEkinhOld = TRUE;
1213 /* for vv, the first half of the integration actually corresponds to the previous step.
1214 So we need information from the last step in the first half of the integration */
1215 if (bGStat || do_per_step(step-1, nstglobalcomm))
1217 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1218 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1219 constr, NULL, FALSE, state->box,
1220 top_global, &bSumEkinhOld,
1223 | (bTemp ? CGLO_TEMPERATURE : 0)
1224 | (bPres ? CGLO_PRESSURE : 0)
1225 | (bPres ? CGLO_CONSTRAINT : 0)
1226 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1227 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1230 /* explanation of above:
1231 a) We compute Ekin at the full time step
1232 if 1) we are using the AveVel Ekin, and it's not the
1233 initial step, or 2) if we are using AveEkin, but need the full
1234 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1235 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1236 EkinAveVel because it's needed for the pressure */
1238 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1243 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1244 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1251 /* We need the kinetic energy at minus the half step for determining
1252 * the full step kinetic energy and possibly for T-coupling.*/
1253 /* This may not be quite working correctly yet . . . . */
1254 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1255 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1256 constr, NULL, FALSE, state->box,
1257 top_global, &bSumEkinhOld,
1258 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1263 if (iterate.bIterationActive &&
1264 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1265 state->veta, &vetanew))
1269 bFirstIterate = FALSE;
1272 if (bTrotter && !bInitStep)
1274 copy_mat(shake_vir, state->svir_prev);
1275 copy_mat(force_vir, state->fvir_prev);
1276 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1278 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1279 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1280 enerd->term[F_EKIN] = trace(ekind->ekin);
1283 /* if it's the initial step, we performed this first step just to get the constraint virial */
1284 if (bInitStep && ir->eI == eiVV)
1286 copy_rvecn(cbuf, state->v, 0, state->natoms);
1290 /* MRS -- now done iterating -- compute the conserved quantity */
1293 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1296 last_ekin = enerd->term[F_EKIN];
1298 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1300 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1302 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1305 sum_dhdl(enerd, state->lambda, ir->fepvals);
1309 /* ######## END FIRST UPDATE STEP ############## */
1310 /* ######## If doing VV, we now have v(dt) ###### */
1313 /* perform extended ensemble sampling in lambda - we don't
1314 actually move to the new state before outputting
1315 statistics, but if performing simulated tempering, we
1316 do update the velocities and the tau_t. */
1318 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, mcrng, state->v, mdatoms);
1319 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1320 copy_df_history(&state_global->dfhist, &state->dfhist);
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, 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 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1343 /* Check whether everything is still allright */
1344 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1345 #ifdef GMX_THREAD_MPI
1350 /* this is just make gs.sig compatible with the hack
1351 of sending signals around by MPI_Reduce with together with
1353 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1355 gs.sig[eglsSTOPCOND] = 1;
1357 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1359 gs.sig[eglsSTOPCOND] = -1;
1361 /* < 0 means stop at next step, > 0 means stop at next NS step */
1365 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1366 gmx_get_signal_name(),
1367 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1371 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1372 gmx_get_signal_name(),
1373 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1375 handled_stop_condition = (int)gmx_get_stop_condition();
1377 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1378 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1379 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1381 /* Signal to terminate the run */
1382 gs.sig[eglsSTOPCOND] = 1;
1385 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1387 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1390 if (bResetCountersHalfMaxH && MASTER(cr) &&
1391 elapsed_time > max_hours*60.0*60.0*0.495)
1393 gs.sig[eglsRESETCOUNTERS] = 1;
1396 if (ir->nstlist == -1 && !bRerunMD)
1398 /* When bGStatEveryStep=FALSE, global_stat is only called
1399 * when we check the atom displacements, not at NS steps.
1400 * This means that also the bonded interaction count check is not
1401 * performed immediately after NS. Therefore a few MD steps could
1402 * be performed with missing interactions.
1403 * But wrong energies are never written to file,
1404 * since energies are only written after global_stat
1407 if (step >= nlh.step_nscheck)
1409 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1410 nlh.scale_tot, state->x);
1414 /* This is not necessarily true,
1415 * but step_nscheck is determined quite conservatively.
1421 /* In parallel we only have to check for checkpointing in steps
1422 * where we do global communication,
1423 * otherwise the other nodes don't know.
1425 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1428 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1429 gs.set[eglsCHKPT] == 0)
1431 gs.sig[eglsCHKPT] = 1;
1434 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1439 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1441 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1443 gmx_bool bIfRandomize;
1444 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1445 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1446 if (constr && bIfRandomize)
1448 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1449 state, fr->bMolPBC, graph, f,
1450 &top->idef, tmp_vir,
1451 cr, nrnb, wcycle, upd, constr,
1452 TRUE, bCalcVir, vetanew);
1457 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1459 gmx_iterate_init(&iterate, TRUE);
1460 /* for iterations, we save these vectors, as we will be redoing the calculations */
1461 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1464 bFirstIterate = TRUE;
1465 while (bFirstIterate || iterate.bIterationActive)
1467 /* We now restore these vectors to redo the calculation with improved extended variables */
1468 if (iterate.bIterationActive)
1470 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1473 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1474 so scroll down for that logic */
1476 /* ######### START SECOND UPDATE STEP ################# */
1477 /* Box is changed in update() when we do pressure coupling,
1478 * but we should still use the old box for energy corrections and when
1479 * writing it to the energy file, so it matches the trajectory files for
1480 * the same timestep above. Make a copy in a separate array.
1482 copy_mat(state->box, lastbox);
1487 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1489 wallcycle_start(wcycle, ewcUPDATE);
1490 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1493 if (iterate.bIterationActive)
1501 /* we use a new value of scalevir to converge the iterations faster */
1502 scalevir = tracevir/trace(shake_vir);
1504 msmul(shake_vir, scalevir, shake_vir);
1505 m_add(force_vir, shake_vir, total_vir);
1506 clear_mat(shake_vir);
1508 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1509 /* We can only do Berendsen coupling after we have summed
1510 * the kinetic energy or virial. Since the happens
1511 * in global_state after update, we should only do it at
1512 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1517 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1518 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1523 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1525 /* velocity half-step update */
1526 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1527 bUpdateDoLR, fr->f_twin, fcd,
1528 ekind, M, upd, FALSE, etrtVELOCITY2,
1529 cr, nrnb, constr, &top->idef);
1532 /* Above, initialize just copies ekinh into ekin,
1533 * it doesn't copy position (for VV),
1534 * and entire integrator for MD.
1537 if (ir->eI == eiVVAK)
1539 copy_rvecn(state->x, cbuf, 0, state->natoms);
1541 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1543 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1544 bUpdateDoLR, fr->f_twin, fcd,
1545 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1546 wallcycle_stop(wcycle, ewcUPDATE);
1548 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1549 fr->bMolPBC, graph, f,
1550 &top->idef, shake_vir,
1551 cr, nrnb, wcycle, upd, constr,
1552 FALSE, bCalcVir, state->veta);
1554 if (ir->eI == eiVVAK)
1556 /* erase F_EKIN and F_TEMP here? */
1557 /* just compute the kinetic energy at the half step to perform a trotter step */
1558 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1559 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1560 constr, NULL, FALSE, lastbox,
1561 top_global, &bSumEkinhOld,
1562 cglo_flags | CGLO_TEMPERATURE
1564 wallcycle_start(wcycle, ewcUPDATE);
1565 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1566 /* now we know the scaling, we can compute the positions again again */
1567 copy_rvecn(cbuf, state->x, 0, state->natoms);
1569 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1571 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1572 bUpdateDoLR, fr->f_twin, fcd,
1573 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1574 wallcycle_stop(wcycle, ewcUPDATE);
1576 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1577 /* are the small terms in the shake_vir here due
1578 * to numerical errors, or are they important
1579 * physically? I'm thinking they are just errors, but not completely sure.
1580 * For now, will call without actually constraining, constr=NULL*/
1581 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1582 state, fr->bMolPBC, graph, f,
1583 &top->idef, tmp_vir,
1584 cr, nrnb, wcycle, upd, NULL,
1590 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1593 if (fr->bSepDVDL && fplog && do_log)
1595 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1599 /* this factor or 2 correction is necessary
1600 because half of the constraint force is removed
1601 in the vv step, so we have to double it. See
1602 the Redmine issue #1255. It is not yet clear
1603 if the factor of 2 is exact, or just a very
1604 good approximation, and this will be
1605 investigated. The next step is to see if this
1606 can be done adding a dhdl contribution from the
1607 rattle step, but this is somewhat more
1608 complicated with the current code. Will be
1609 investigated, hopefully for 4.6.3. However,
1610 this current solution is much better than
1611 having it completely wrong.
1613 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1617 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1622 /* Need to unshift here */
1623 unshift_self(graph, state->box, state->x);
1628 wallcycle_start(wcycle, ewcVSITECONSTR);
1631 shift_self(graph, state->box, state->x);
1633 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1634 top->idef.iparams, top->idef.il,
1635 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1639 unshift_self(graph, state->box, state->x);
1641 wallcycle_stop(wcycle, ewcVSITECONSTR);
1644 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1645 /* With Leap-Frog we can skip compute_globals at
1646 * non-communication steps, but we need to calculate
1647 * the kinetic energy one step before communication.
1649 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1651 if (ir->nstlist == -1 && bFirstIterate)
1653 gs.sig[eglsNABNSB] = nlh.nabnsb;
1655 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1656 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1658 bFirstIterate ? &gs : NULL,
1659 (step_rel % gs.nstms == 0) &&
1660 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1662 top_global, &bSumEkinhOld,
1664 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1665 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1666 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1667 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1668 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1669 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1672 if (ir->nstlist == -1 && bFirstIterate)
1674 nlh.nabnsb = gs.set[eglsNABNSB];
1675 gs.set[eglsNABNSB] = 0;
1678 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1679 /* ############# END CALC EKIN AND PRESSURE ################# */
1681 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1682 the virial that should probably be addressed eventually. state->veta has better properies,
1683 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1684 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1686 if (iterate.bIterationActive &&
1687 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1688 trace(shake_vir), &tracevir))
1692 bFirstIterate = FALSE;
1695 if (!bVV || bRerunMD)
1697 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1698 sum_dhdl(enerd, state->lambda, ir->fepvals);
1700 update_box(fplog, step, ir, mdatoms, state, f,
1701 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1703 /* ################# END UPDATE STEP 2 ################# */
1704 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1706 /* The coordinates (x) were unshifted in update */
1709 /* We will not sum ekinh_old,
1710 * so signal that we still have to do it.
1712 bSumEkinhOld = TRUE;
1715 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1717 /* use the directly determined last velocity, not actually the averaged half steps */
1718 if (bTrotter && ir->eI == eiVV)
1720 enerd->term[F_EKIN] = last_ekin;
1722 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1726 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1730 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1732 /* ######### END PREPARING EDR OUTPUT ########### */
1737 gmx_bool do_dr, do_or;
1739 if (fplog && do_log && bDoExpanded)
1741 /* only needed if doing expanded ensemble */
1742 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1743 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1745 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1749 upd_mdebin(mdebin, bDoDHDL, TRUE,
1750 t, mdatoms->tmass, enerd, state,
1751 ir->fepvals, ir->expandedvals, lastbox,
1752 shake_vir, force_vir, total_vir, pres,
1753 ekind, mu_tot, constr);
1757 upd_mdebin_step(mdebin);
1760 do_dr = do_per_step(step, ir->nstdisreout);
1761 do_or = do_per_step(step, ir->nstorireout);
1763 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
1765 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1767 if (ir->ePull != epullNO)
1769 pull_print_output(ir->pull, step, t);
1772 if (do_per_step(step, ir->nstlog))
1774 if (fflush(fplog) != 0)
1776 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1782 /* Have to do this part _after_ outputting the logfile and the edr file */
1783 /* Gets written into the state at the beginning of next loop*/
1784 state->fep_state = lamnew;
1786 /* Print the remaining wall clock time for the run */
1787 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1791 fprintf(stderr, "\n");
1793 print_time(stderr, walltime_accounting, step, ir, cr);
1796 /* Replica exchange */
1798 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1799 do_per_step(step, repl_ex_nst))
1801 bExchanged = replica_exchange(fplog, cr, repl_ex,
1802 state_global, enerd,
1805 if (bExchanged && DOMAINDECOMP(cr))
1807 dd_partition_system(fplog, step, cr, TRUE, 1,
1808 state_global, top_global, ir,
1809 state, &f, mdatoms, top, fr,
1810 vsite, shellfc, constr,
1811 nrnb, wcycle, FALSE);
1817 bStartingFromCpt = FALSE;
1819 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1820 /* With all integrators, except VV, we need to retain the pressure
1821 * at the current step for coupling at the next step.
1823 if ((state->flags & (1<<estPRES_PREV)) &&
1825 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1827 /* Store the pressure in t_state for pressure coupling
1828 * at the next MD step.
1830 copy_mat(pres, state->pres_prev);
1833 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1835 if ( (membed != NULL) && (!bLastStep) )
1837 rescale_membed(step_rel, membed, state_global->x);
1844 /* read next frame from input trajectory */
1845 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1850 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1854 if (!bRerunMD || !rerun_fr.bStep)
1856 /* increase the MD step number */
1861 cycles = wallcycle_stop(wcycle, ewcSTEP);
1862 if (DOMAINDECOMP(cr) && wcycle)
1864 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1867 if (bPMETuneRunning || bPMETuneTry)
1869 /* PME grid + cut-off optimization with GPUs or PME nodes */
1871 /* Count the total cycles over the last steps */
1872 cycles_pmes += cycles;
1874 /* We can only switch cut-off at NS steps */
1875 if (step % ir->nstlist == 0)
1877 /* PME grid + cut-off optimization with GPUs or PME nodes */
1880 if (DDMASTER(cr->dd))
1882 /* PME node load is too high, start tuning */
1883 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1885 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1887 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1889 bPMETuneTry = FALSE;
1892 if (bPMETuneRunning)
1894 /* init_step might not be a multiple of nstlist,
1895 * but the first cycle is always skipped anyhow.
1898 pme_load_balance(pme_loadbal, cr,
1899 (bVerbose && MASTER(cr)) ? stderr : NULL,
1901 ir, state, cycles_pmes,
1902 fr->ic, fr->nbv, &fr->pmedata,
1905 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1906 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1907 fr->rlist = fr->ic->rlist;
1908 fr->rlistlong = fr->ic->rlistlong;
1909 fr->rcoulomb = fr->ic->rcoulomb;
1910 fr->rvdw = fr->ic->rvdw;
1916 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1917 gs.set[eglsRESETCOUNTERS] != 0)
1919 /* Reset all the counters related to performance over the run */
1920 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1921 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1922 wcycle_set_reset_counters(wcycle, -1);
1923 if (!(cr->duty & DUTY_PME))
1925 /* Tell our PME node to reset its counters */
1926 gmx_pme_send_resetcounters(cr, step);
1928 /* Correct max_hours for the elapsed time */
1929 max_hours -= elapsed_time/(60.0*60.0);
1930 bResetCountersHalfMaxH = FALSE;
1931 gs.set[eglsRESETCOUNTERS] = 0;
1935 /* End of main MD loop */
1938 /* Stop measuring walltime */
1939 walltime_accounting_end(walltime_accounting);
1941 if (bRerunMD && MASTER(cr))
1946 if (!(cr->duty & DUTY_PME))
1948 /* Tell the PME only node to finish */
1949 gmx_pme_send_finish(cr);
1954 if (ir->nstcalcenergy > 0 && !bRerunMD)
1956 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
1957 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1965 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1967 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)));
1968 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1971 if (pme_loadbal != NULL)
1973 pme_loadbal_done(pme_loadbal, cr, fplog,
1974 fr->nbv != NULL && fr->nbv->bUseGPU);
1977 if (shellfc && fplog)
1979 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1980 (nconverged*100.0)/step_rel);
1981 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1985 if (repl_ex_nst > 0 && MASTER(cr))
1987 print_replica_exchange_statistics(fplog, repl_ex);
1990 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);