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.
55 #include "md_support.h"
56 #include "md_logging.h"
70 #include "domdec_network.h"
76 #include "compute_io.h"
78 #include "checkpoint.h"
79 #include "mtop_util.h"
80 #include "sighandler.h"
83 #include "pme_loadbal.h"
86 #include "types/nlistheuristics.h"
87 #include "types/iteratedconstraints.h"
88 #include "nbnxn_cuda_data_mgmt.h"
90 #include "gromacs/utility/gmxmpi.h"
91 #include "gromacs/fileio/confio.h"
92 #include "gromacs/fileio/trajectory_writing.h"
93 #include "gromacs/fileio/trnio.h"
94 #include "gromacs/fileio/trxio.h"
95 #include "gromacs/fileio/xtcio.h"
96 #include "gromacs/timing/wallcycle.h"
97 #include "gromacs/timing/walltime_accounting.h"
100 #include "corewrap.h"
103 static void reset_all_counters(FILE *fplog, t_commrec *cr,
104 gmx_large_int_t step,
105 gmx_large_int_t *step_rel, t_inputrec *ir,
106 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
107 gmx_walltime_accounting_t walltime_accounting,
108 nbnxn_cuda_ptr_t cu_nbv)
110 char sbuf[STEPSTRSIZE];
112 /* Reset all the counters related to performance over the run */
113 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
114 gmx_step_str(step, sbuf));
118 nbnxn_cuda_reset_timings(cu_nbv);
121 wallcycle_stop(wcycle, ewcRUN);
122 wallcycle_reset_all(wcycle);
123 if (DOMAINDECOMP(cr))
125 reset_dd_statistics_counters(cr->dd);
128 ir->init_step += *step_rel;
129 ir->nsteps -= *step_rel;
131 wallcycle_start(wcycle, ewcRUN);
132 walltime_accounting_start(walltime_accounting);
133 print_date_and_time(fplog, cr->nodeid, "Restarted time", walltime_accounting);
136 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
137 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
139 gmx_vsite_t *vsite, gmx_constr_t constr,
140 int stepout, t_inputrec *ir,
141 gmx_mtop_t *top_global,
143 t_state *state_global,
145 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
146 gmx_edsam_t ed, t_forcerec *fr,
147 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
148 real cpt_period, real max_hours,
149 const char gmx_unused *deviceOptions,
151 gmx_walltime_accounting_t walltime_accounting)
154 gmx_large_int_t step, step_rel;
156 double t, t0, lam0[efptNR];
157 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
158 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
159 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
160 bBornRadii, bStartingFromCpt;
161 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
162 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
163 bForceUpdate = FALSE, bCPT;
164 gmx_bool bMasterState;
165 int force_flags, cglo_flags;
166 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
171 t_state *bufstate = NULL;
172 matrix *scale_tot, pcoupl_mu, M, ebox;
175 gmx_repl_ex_t repl_ex = NULL;
178 t_mdebin *mdebin = NULL;
179 t_state *state = NULL;
180 rvec *f_global = NULL;
181 gmx_enerdata_t *enerd;
183 gmx_global_stat_t gstat;
184 gmx_update_t upd = NULL;
185 t_graph *graph = NULL;
187 gmx_rng_t mcrng = NULL;
188 gmx_groups_t *groups;
189 gmx_ekindata_t *ekind, *ekind_save;
190 gmx_shellfc_t shellfc;
191 int count, nconverged = 0;
194 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
196 gmx_bool bResetCountersHalfMaxH = FALSE;
197 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
198 gmx_bool bUpdateDoLR;
203 real veta_save, scalevir, tracevir;
209 real saved_conserved_quantity = 0;
214 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
215 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
216 gmx_iterate_t iterate;
217 gmx_large_int_t multisim_nsteps = -1; /* number of steps to do before first multisim
218 simulation stops. If equal to zero, don't
219 communicate any more between multisims.*/
220 /* PME load balancing data for GPU kernels */
221 pme_load_balancing_t pme_loadbal = NULL;
223 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
226 /* Temporary addition for FAHCORE checkpointing */
230 /* Check for special mdrun options */
231 bRerunMD = (Flags & MD_RERUN);
232 bAppend = (Flags & MD_APPENDFILES);
233 if (Flags & MD_RESETCOUNTERSHALFWAY)
237 /* Signal to reset the counters half the simulation steps. */
238 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
240 /* Signal to reset the counters halfway the simulation time. */
241 bResetCountersHalfMaxH = (max_hours > 0);
244 /* md-vv uses averaged full step velocities for T-control
245 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
246 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
248 if (bVV) /* to store the initial velocities while computing virial */
250 snew(cbuf, top_global->natoms);
252 /* all the iteratative cases - only if there are constraints */
253 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
254 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
255 false in this step. The correct value, true or false,
256 is set at each step, as it depends on the frequency of temperature
257 and pressure control.*/
258 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
262 /* Since we don't know if the frames read are related in any way,
263 * rebuild the neighborlist at every step.
266 ir->nstcalcenergy = 1;
270 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
272 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
273 bGStatEveryStep = (nstglobalcomm == 1);
275 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
278 "To reduce the energy communication with nstlist = -1\n"
279 "the neighbor list validity should not be checked at every step,\n"
280 "this means that exact integration is not guaranteed.\n"
281 "The neighbor list validity is checked after:\n"
282 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
283 "In most cases this will result in exact integration.\n"
284 "This reduces the energy communication by a factor of 2 to 3.\n"
285 "If you want less energy communication, set nstlist > 3.\n\n");
292 groups = &top_global->groups;
295 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
296 &(state_global->fep_state), lam0,
297 nrnb, top_global, &upd,
298 nfile, fnm, &outf, &mdebin,
299 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
301 clear_mat(total_vir);
303 /* Energy terms and groups */
305 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
307 if (DOMAINDECOMP(cr))
313 snew(f, top_global->natoms);
316 /* Kinetic energy data */
318 init_ekindata(fplog, top_global, &(ir->opts), ekind);
319 /* needed for iteration of constraints */
321 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
322 /* Copy the cos acceleration to the groups struct */
323 ekind->cosacc.cos_accel = ir->cos_accel;
325 gstat = global_stat_init(ir);
328 /* Check for polarizable models and flexible constraints */
329 shellfc = init_shell_flexcon(fplog,
330 top_global, n_flexible_constraints(constr),
331 (ir->bContinuation ||
332 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
333 NULL : state_global->x);
337 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 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
345 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
346 if ((io > 2000) && MASTER(cr))
349 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
354 if (DOMAINDECOMP(cr))
356 top = dd_init_local_top(top_global);
359 dd_init_local_state(cr->dd, state_global, state);
361 if (DDMASTER(cr->dd) && ir->nstfout)
363 snew(f_global, state_global->natoms);
370 /* Initialize the particle decomposition and split the topology */
371 top = split_system(fplog, top_global, ir, cr);
373 pd_cg_range(cr, &fr->cg0, &fr->hcg);
374 pd_at_range(cr, &a0, &a1);
378 top = gmx_mtop_generate_local_top(top_global, ir);
381 a1 = top_global->natoms;
384 forcerec_set_excl_load(fr, top, cr);
386 state = partdec_init_local_state(cr, state_global);
389 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
393 set_vsite_top(vsite, top, mdatoms, cr);
396 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
398 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
403 make_local_shells(cr, mdatoms, shellfc);
406 setup_bonded_threading(fr, &top->idef);
408 if (ir->pull && PAR(cr))
410 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
414 if (DOMAINDECOMP(cr))
416 /* Distribute the charge groups over the nodes from the master node */
417 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
418 state_global, top_global, ir,
419 state, &f, mdatoms, top, fr,
420 vsite, shellfc, constr,
421 nrnb, wcycle, FALSE);
425 update_mdatoms(mdatoms, state->lambda[efptMASS]);
427 if (opt2bSet("-cpi", nfile, fnm))
429 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
433 bStateFromCP = FALSE;
438 init_expanded_ensemble(bStateFromCP,ir,&mcrng,&state->dfhist);
445 /* Update mdebin with energy history if appending to output files */
446 if (Flags & MD_APPENDFILES)
448 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
452 /* We might have read an energy history from checkpoint,
453 * free the allocated memory and reset the counts.
455 done_energyhistory(&state_global->enerhist);
456 init_energyhistory(&state_global->enerhist);
459 /* Set the initial energy history in state by updating once */
460 update_energyhistory(&state_global->enerhist, mdebin);
463 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
465 /* Set the random state if we read a checkpoint file */
466 set_stochd_state(upd, state);
469 if (state->flags & (1<<estMC_RNG))
471 set_mc_state(mcrng, state);
474 /* Initialize constraints */
477 if (!DOMAINDECOMP(cr))
479 set_constraints(constr, top, ir, mdatoms, cr);
485 /* We need to be sure replica exchange can only occur
486 * when the energies are current */
487 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
488 "repl_ex_nst", &repl_ex_nst);
489 /* This check needs to happen before inter-simulation
490 * signals are initialized, too */
492 if (repl_ex_nst > 0 && MASTER(cr))
494 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
495 repl_ex_nst, repl_ex_nex, repl_ex_seed);
498 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
499 * With perturbed charges with soft-core we should not change the cut-off.
501 if ((Flags & MD_TUNEPME) &&
502 EEL_PME(fr->eeltype) &&
503 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
504 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
507 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
509 if (cr->duty & DUTY_PME)
511 /* Start tuning right away, as we can't measure the load */
512 bPMETuneRunning = TRUE;
516 /* Separate PME nodes, we can measure the PP/PME load balance */
521 if (!ir->bContinuation && !bRerunMD)
523 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
525 /* Set the velocities of frozen particles to zero */
526 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
528 for (m = 0; m < DIM; m++)
530 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
540 /* Constrain the initial coordinates and velocities */
541 do_constrain_first(fplog, constr, ir, mdatoms, state,
546 /* Construct the virtual sites for the initial configuration */
547 construct_vsites(vsite, state->x, ir->delta_t, NULL,
548 top->idef.iparams, top->idef.il,
549 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
555 /* set free energy calculation frequency as the minimum
556 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
557 nstfep = ir->fepvals->nstdhdl;
560 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
564 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
567 /* I'm assuming we need global communication the first time! MRS */
568 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
569 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
570 | (bVV ? CGLO_PRESSURE : 0)
571 | (bVV ? CGLO_CONSTRAINT : 0)
572 | (bRerunMD ? CGLO_RERUNMD : 0)
573 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
575 bSumEkinhOld = FALSE;
576 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
577 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
578 constr, NULL, FALSE, state->box,
579 top_global, &bSumEkinhOld, cglo_flags);
580 if (ir->eI == eiVVAK)
582 /* a second call to get the half step temperature initialized as well */
583 /* we do the same call as above, but turn the pressure off -- internally to
584 compute_globals, this is recognized as a velocity verlet half-step
585 kinetic energy calculation. This minimized excess variables, but
586 perhaps loses some logic?*/
588 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
589 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
590 constr, NULL, FALSE, state->box,
591 top_global, &bSumEkinhOld,
592 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
595 /* Calculate the initial half step temperature, and save the ekinh_old */
596 if (!(Flags & MD_STARTFROMCPT))
598 for (i = 0; (i < ir->opts.ngtc); i++)
600 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
605 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
606 and there is no previous step */
609 /* if using an iterative algorithm, we need to create a working directory for the state. */
612 bufstate = init_bufstate(state);
615 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
616 temperature control */
617 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
621 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
624 "RMS relative constraint deviation after constraining: %.2e\n",
625 constr_rmsd(constr, FALSE));
627 if (EI_STATE_VELOCITY(ir->eI))
629 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
633 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
634 " input trajectory '%s'\n\n",
635 *(top_global->name), opt2fn("-rerun", nfile, fnm));
638 fprintf(stderr, "Calculated time to finish depends on nsteps from "
639 "run input file,\nwhich may not correspond to the time "
640 "needed to process input trajectory.\n\n");
646 fprintf(stderr, "starting mdrun '%s'\n",
647 *(top_global->name));
650 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
654 sprintf(tbuf, "%s", "infinite");
656 if (ir->init_step > 0)
658 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
659 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
660 gmx_step_str(ir->init_step, sbuf2),
661 ir->init_step*ir->delta_t);
665 fprintf(stderr, "%s steps, %s ps.\n",
666 gmx_step_str(ir->nsteps, sbuf), tbuf);
669 fprintf(fplog, "\n");
672 /* Set and write start time */
673 walltime_accounting_start(walltime_accounting);
674 print_date_and_time(fplog, cr->nodeid, "Started mdrun", walltime_accounting);
675 wallcycle_start(wcycle, ewcRUN);
678 fprintf(fplog, "\n");
681 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
683 chkpt_ret = fcCheckPointParallel( cr->nodeid,
687 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
692 /***********************************************************
696 ************************************************************/
698 /* if rerunMD then read coordinates and velocities from input trajectory */
701 if (getenv("GMX_FORCE_UPDATE"))
709 bNotLastFrame = read_first_frame(oenv, &status,
710 opt2fn("-rerun", nfile, fnm),
711 &rerun_fr, TRX_NEED_X | TRX_READ_V);
712 if (rerun_fr.natoms != top_global->natoms)
715 "Number of atoms in trajectory (%d) does not match the "
716 "run input file (%d)\n",
717 rerun_fr.natoms, top_global->natoms);
719 if (ir->ePBC != epbcNONE)
723 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);
725 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
727 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
734 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
737 if (ir->ePBC != epbcNONE)
739 /* Set the shift vectors.
740 * Necessary here when have a static box different from the tpr box.
742 calc_shifts(rerun_fr.box, fr->shift_vec);
746 /* loop over MD steps or if rerunMD to end of input trajectory */
748 /* Skip the first Nose-Hoover integration when we get the state from tpx */
749 bStateFromTPX = !bStateFromCP;
750 bInitStep = bFirstStep && (bStateFromTPX || bVV);
751 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
753 bSumEkinhOld = FALSE;
756 init_global_signals(&gs, cr, ir, repl_ex_nst);
758 step = ir->init_step;
761 if (ir->nstlist == -1)
763 init_nlistheuristics(&nlh, bGStatEveryStep, step);
766 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
768 /* check how many steps are left in other sims */
769 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
773 /* and stop now if we should */
774 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
775 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
776 while (!bLastStep || (bRerunMD && bNotLastFrame))
779 wallcycle_start(wcycle, ewcSTEP);
785 step = rerun_fr.step;
786 step_rel = step - ir->init_step;
799 bLastStep = (step_rel == ir->nsteps);
800 t = t0 + step*ir->delta_t;
803 if (ir->efep != efepNO || ir->bSimTemp)
805 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
806 requiring different logic. */
808 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
809 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
810 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
811 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
812 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
817 update_annealing_target_temp(&(ir->opts), t);
822 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
824 for (i = 0; i < state_global->natoms; i++)
826 copy_rvec(rerun_fr.x[i], state_global->x[i]);
830 for (i = 0; i < state_global->natoms; i++)
832 copy_rvec(rerun_fr.v[i], state_global->v[i]);
837 for (i = 0; i < state_global->natoms; i++)
839 clear_rvec(state_global->v[i]);
843 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
844 " Ekin, temperature and pressure are incorrect,\n"
845 " the virial will be incorrect when constraints are present.\n"
847 bRerunWarnNoV = FALSE;
851 copy_mat(rerun_fr.box, state_global->box);
852 copy_mat(state_global->box, state->box);
854 if (vsite && (Flags & MD_RERUN_VSITE))
856 if (DOMAINDECOMP(cr))
858 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
862 /* Following is necessary because the graph may get out of sync
863 * with the coordinates if we only have every N'th coordinate set
865 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
866 shift_self(graph, state->box, state->x);
868 construct_vsites(vsite, state->x, ir->delta_t, state->v,
869 top->idef.iparams, top->idef.il,
870 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
873 unshift_self(graph, state->box, state->x);
878 /* Stop Center of Mass motion */
879 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
883 /* for rerun MD always do Neighbour Searching */
884 bNS = (bFirstStep || ir->nstlist != 0);
889 /* Determine whether or not to do Neighbour Searching and LR */
890 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
892 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
893 (ir->nstlist == -1 && nlh.nabnsb > 0));
895 if (bNS && ir->nstlist == -1)
897 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
901 /* check whether we should stop because another simulation has
905 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
906 (multisim_nsteps != ir->nsteps) )
913 "Stopping simulation %d because another one has finished\n",
917 gs.sig[eglsCHKPT] = 1;
922 /* < 0 means stop at next step, > 0 means stop at next NS step */
923 if ( (gs.set[eglsSTOPCOND] < 0) ||
924 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
929 /* Determine whether or not to update the Born radii if doing GB */
930 bBornRadii = bFirstStep;
931 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
936 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
937 do_verbose = bVerbose &&
938 (step % stepout == 0 || bFirstStep || bLastStep);
940 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
948 bMasterState = FALSE;
949 /* Correct the new box if it is too skewed */
950 if (DYNAMIC_BOX(*ir))
952 if (correct_box(fplog, step, state->box, graph))
957 if (DOMAINDECOMP(cr) && bMasterState)
959 dd_collect_state(cr->dd, state, state_global);
963 if (DOMAINDECOMP(cr))
965 /* Repartition the domain decomposition */
966 wallcycle_start(wcycle, ewcDOMDEC);
967 dd_partition_system(fplog, step, cr,
968 bMasterState, nstglobalcomm,
969 state_global, top_global, ir,
970 state, &f, mdatoms, top, fr,
971 vsite, shellfc, constr,
973 do_verbose && !bPMETuneRunning);
974 wallcycle_stop(wcycle, ewcDOMDEC);
975 /* If using an iterative integrator, reallocate space to match the decomposition */
979 if (MASTER(cr) && do_log)
981 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
984 if (ir->efep != efepNO)
986 update_mdatoms(mdatoms, state->lambda[efptMASS]);
989 if ((bRerunMD && rerun_fr.bV) || bExchanged)
992 /* We need the kinetic energy at minus the half step for determining
993 * the full step kinetic energy and possibly for T-coupling.*/
994 /* This may not be quite working correctly yet . . . . */
995 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
996 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
997 constr, NULL, FALSE, state->box,
998 top_global, &bSumEkinhOld,
999 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1001 clear_mat(force_vir);
1003 /* We write a checkpoint at this MD step when:
1004 * either at an NS step when we signalled through gs,
1005 * or at the last step (but not when we do not want confout),
1006 * but never at the first step or with rerun.
1008 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1009 (bLastStep && (Flags & MD_CONFOUT))) &&
1010 step > ir->init_step && !bRerunMD);
1013 gs.set[eglsCHKPT] = 0;
1016 /* Determine the energy and pressure:
1017 * at nstcalcenergy steps and at energy output steps (set below).
1019 if (EI_VV(ir->eI) && (!bInitStep))
1021 /* for vv, the first half of the integration actually corresponds
1022 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1023 but the virial needs to be calculated on both the current step and the 'next' step. Future
1024 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1026 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1027 bCalcVir = bCalcEner ||
1028 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1032 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1033 bCalcVir = bCalcEner ||
1034 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1037 /* Do we need global communication ? */
1038 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1039 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1040 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1042 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1044 if (do_ene || do_log)
1051 /* these CGLO_ options remain the same throughout the iteration */
1052 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1053 (bGStat ? CGLO_GSTAT : 0)
1056 force_flags = (GMX_FORCE_STATECHANGED |
1057 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1058 GMX_FORCE_ALLFORCES |
1060 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1061 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1062 (bDoFEP ? GMX_FORCE_DHDL : 0)
1067 if (do_per_step(step, ir->nstcalclr))
1069 force_flags |= GMX_FORCE_DO_LR;
1075 /* Now is the time to relax the shells */
1076 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1077 ir, bNS, force_flags,
1080 state, f, force_vir, mdatoms,
1081 nrnb, wcycle, graph, groups,
1082 shellfc, fr, bBornRadii, t, mu_tot,
1094 /* The coordinates (x) are shifted (to get whole molecules)
1096 * This is parallellized as well, and does communication too.
1097 * Check comments in sim_util.c
1099 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1100 state->box, state->x, &state->hist,
1101 f, force_vir, mdatoms, enerd, fcd,
1102 state->lambda, graph,
1103 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1104 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1107 if (bVV && !bStartingFromCpt && !bRerunMD)
1108 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1110 if (ir->eI == eiVV && bInitStep)
1112 /* if using velocity verlet with full time step Ekin,
1113 * take the first half step only to compute the
1114 * virial for the first step. From there,
1115 * revert back to the initial coordinates
1116 * so that the input is actually the initial step.
1118 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1122 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1123 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1126 /* If we are using twin-range interactions where the long-range component
1127 * is only evaluated every nstcalclr>1 steps, we should do a special update
1128 * step to combine the long-range forces on these steps.
1129 * For nstcalclr=1 this is not done, since the forces would have been added
1130 * directly to the short-range forces already.
1132 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1134 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1135 f, bUpdateDoLR, fr->f_twin, fcd,
1136 ekind, M, upd, bInitStep, etrtVELOCITY1,
1137 cr, nrnb, constr, &top->idef);
1139 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1141 gmx_iterate_init(&iterate, TRUE);
1143 /* for iterations, we save these vectors, as we will be self-consistently iterating
1146 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1148 /* save the state */
1149 if (iterate.bIterationActive)
1151 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1154 bFirstIterate = TRUE;
1155 while (bFirstIterate || iterate.bIterationActive)
1157 if (iterate.bIterationActive)
1159 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1160 if (bFirstIterate && bTrotter)
1162 /* The first time through, we need a decent first estimate
1163 of veta(t+dt) to compute the constraints. Do
1164 this by computing the box volume part of the
1165 trotter integration at this time. Nothing else
1166 should be changed by this routine here. If
1167 !(first time), we start with the previous value
1170 veta_save = state->veta;
1171 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1172 vetanew = state->veta;
1173 state->veta = veta_save;
1178 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1180 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1181 state, fr->bMolPBC, graph, f,
1182 &top->idef, shake_vir,
1183 cr, nrnb, wcycle, upd, constr,
1184 TRUE, bCalcVir, vetanew);
1188 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1194 /* Need to unshift here if a do_force has been
1195 called in the previous step */
1196 unshift_self(graph, state->box, state->x);
1199 /* if VV, compute the pressure and constraints */
1200 /* For VV2, we strictly only need this if using pressure
1201 * control, but we really would like to have accurate pressures
1203 * Think about ways around this in the future?
1204 * For now, keep this choice in comments.
1206 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1207 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1209 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1210 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1212 bSumEkinhOld = TRUE;
1214 /* for vv, the first half of the integration actually corresponds to the previous step.
1215 So we need information from the last step in the first half of the integration */
1216 if (bGStat || do_per_step(step-1, nstglobalcomm))
1218 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1219 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1220 constr, NULL, FALSE, state->box,
1221 top_global, &bSumEkinhOld,
1224 | (bTemp ? CGLO_TEMPERATURE : 0)
1225 | (bPres ? CGLO_PRESSURE : 0)
1226 | (bPres ? CGLO_CONSTRAINT : 0)
1227 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1228 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1231 /* explanation of above:
1232 a) We compute Ekin at the full time step
1233 if 1) we are using the AveVel Ekin, and it's not the
1234 initial step, or 2) if we are using AveEkin, but need the full
1235 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1236 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1237 EkinAveVel because it's needed for the pressure */
1239 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1244 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1245 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1252 /* We need the kinetic energy at minus the half step for determining
1253 * the full step kinetic energy and possibly for T-coupling.*/
1254 /* This may not be quite working correctly yet . . . . */
1255 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1256 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1257 constr, NULL, FALSE, state->box,
1258 top_global, &bSumEkinhOld,
1259 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1264 if (iterate.bIterationActive &&
1265 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1266 state->veta, &vetanew))
1270 bFirstIterate = FALSE;
1273 if (bTrotter && !bInitStep)
1275 copy_mat(shake_vir, state->svir_prev);
1276 copy_mat(force_vir, state->fvir_prev);
1277 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1279 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1280 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1281 enerd->term[F_EKIN] = trace(ekind->ekin);
1284 /* if it's the initial step, we performed this first step just to get the constraint virial */
1285 if (bInitStep && ir->eI == eiVV)
1287 copy_rvecn(cbuf, state->v, 0, state->natoms);
1291 /* MRS -- now done iterating -- compute the conserved quantity */
1294 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1297 last_ekin = enerd->term[F_EKIN];
1299 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1301 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1303 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1306 sum_dhdl(enerd, state->lambda, ir->fepvals);
1310 /* ######## END FIRST UPDATE STEP ############## */
1311 /* ######## If doing VV, we now have v(dt) ###### */
1314 /* perform extended ensemble sampling in lambda - we don't
1315 actually move to the new state before outputting
1316 statistics, but if performing simulated tempering, we
1317 do update the velocities and the tau_t. */
1319 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, mcrng, state->v, mdatoms);
1320 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1321 copy_df_history(&state_global->dfhist,&state->dfhist);
1324 /* Now we have the energies and forces corresponding to the
1325 * coordinates at time t. We must output all of this before
1328 do_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1329 ir, state, state_global, top_global, fr, upd,
1330 outf, mdebin, ekind, f, f_global,
1331 wcycle, mcrng, &nchkpt,
1332 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1335 /* kludge -- virial is lost with restart for NPT control. Must restart */
1336 if (bStartingFromCpt && bVV)
1338 copy_mat(state->svir_prev, shake_vir);
1339 copy_mat(state->fvir_prev, force_vir);
1342 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1344 /* Check whether everything is still allright */
1345 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1346 #ifdef GMX_THREAD_MPI
1351 /* this is just make gs.sig compatible with the hack
1352 of sending signals around by MPI_Reduce with together with
1354 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1356 gs.sig[eglsSTOPCOND] = 1;
1358 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1360 gs.sig[eglsSTOPCOND] = -1;
1362 /* < 0 means stop at next step, > 0 means stop at next NS step */
1366 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1367 gmx_get_signal_name(),
1368 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1372 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1373 gmx_get_signal_name(),
1374 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1376 handled_stop_condition = (int)gmx_get_stop_condition();
1378 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1379 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1380 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1382 /* Signal to terminate the run */
1383 gs.sig[eglsSTOPCOND] = 1;
1386 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1388 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1391 if (bResetCountersHalfMaxH && MASTER(cr) &&
1392 elapsed_time > max_hours*60.0*60.0*0.495)
1394 gs.sig[eglsRESETCOUNTERS] = 1;
1397 if (ir->nstlist == -1 && !bRerunMD)
1399 /* When bGStatEveryStep=FALSE, global_stat is only called
1400 * when we check the atom displacements, not at NS steps.
1401 * This means that also the bonded interaction count check is not
1402 * performed immediately after NS. Therefore a few MD steps could
1403 * be performed with missing interactions.
1404 * But wrong energies are never written to file,
1405 * since energies are only written after global_stat
1408 if (step >= nlh.step_nscheck)
1410 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1411 nlh.scale_tot, state->x);
1415 /* This is not necessarily true,
1416 * but step_nscheck is determined quite conservatively.
1422 /* In parallel we only have to check for checkpointing in steps
1423 * where we do global communication,
1424 * otherwise the other nodes don't know.
1426 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1429 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1430 gs.set[eglsCHKPT] == 0)
1432 gs.sig[eglsCHKPT] = 1;
1435 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1440 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1442 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1444 gmx_bool bIfRandomize;
1445 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1446 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1447 if (constr && bIfRandomize)
1449 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1450 state, fr->bMolPBC, graph, f,
1451 &top->idef, tmp_vir,
1452 cr, nrnb, wcycle, upd, constr,
1453 TRUE, bCalcVir, vetanew);
1458 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1460 gmx_iterate_init(&iterate, TRUE);
1461 /* for iterations, we save these vectors, as we will be redoing the calculations */
1462 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1465 bFirstIterate = TRUE;
1466 while (bFirstIterate || iterate.bIterationActive)
1468 /* We now restore these vectors to redo the calculation with improved extended variables */
1469 if (iterate.bIterationActive)
1471 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1474 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1475 so scroll down for that logic */
1477 /* ######### START SECOND UPDATE STEP ################# */
1478 /* Box is changed in update() when we do pressure coupling,
1479 * but we should still use the old box for energy corrections and when
1480 * writing it to the energy file, so it matches the trajectory files for
1481 * the same timestep above. Make a copy in a separate array.
1483 copy_mat(state->box, lastbox);
1488 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1490 wallcycle_start(wcycle, ewcUPDATE);
1491 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1494 if (iterate.bIterationActive)
1502 /* we use a new value of scalevir to converge the iterations faster */
1503 scalevir = tracevir/trace(shake_vir);
1505 msmul(shake_vir, scalevir, shake_vir);
1506 m_add(force_vir, shake_vir, total_vir);
1507 clear_mat(shake_vir);
1509 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1510 /* We can only do Berendsen coupling after we have summed
1511 * the kinetic energy or virial. Since the happens
1512 * in global_state after update, we should only do it at
1513 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1518 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1519 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1524 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1526 /* velocity half-step update */
1527 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1528 bUpdateDoLR, fr->f_twin, fcd,
1529 ekind, M, upd, FALSE, etrtVELOCITY2,
1530 cr, nrnb, constr, &top->idef);
1533 /* Above, initialize just copies ekinh into ekin,
1534 * it doesn't copy position (for VV),
1535 * and entire integrator for MD.
1538 if (ir->eI == eiVVAK)
1540 copy_rvecn(state->x, cbuf, 0, state->natoms);
1542 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1544 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1545 bUpdateDoLR, fr->f_twin, fcd,
1546 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1547 wallcycle_stop(wcycle, ewcUPDATE);
1549 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1550 fr->bMolPBC, graph, f,
1551 &top->idef, shake_vir,
1552 cr, nrnb, wcycle, upd, constr,
1553 FALSE, bCalcVir, state->veta);
1555 if (ir->eI == eiVVAK)
1557 /* erase F_EKIN and F_TEMP here? */
1558 /* just compute the kinetic energy at the half step to perform a trotter step */
1559 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1560 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1561 constr, NULL, FALSE, lastbox,
1562 top_global, &bSumEkinhOld,
1563 cglo_flags | CGLO_TEMPERATURE
1565 wallcycle_start(wcycle, ewcUPDATE);
1566 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1567 /* now we know the scaling, we can compute the positions again again */
1568 copy_rvecn(cbuf, state->x, 0, state->natoms);
1570 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1572 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1573 bUpdateDoLR, fr->f_twin, fcd,
1574 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1575 wallcycle_stop(wcycle, ewcUPDATE);
1577 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1578 /* are the small terms in the shake_vir here due
1579 * to numerical errors, or are they important
1580 * physically? I'm thinking they are just errors, but not completely sure.
1581 * For now, will call without actually constraining, constr=NULL*/
1582 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1583 state, fr->bMolPBC, graph, f,
1584 &top->idef, tmp_vir,
1585 cr, nrnb, wcycle, upd, NULL,
1591 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1594 if (fr->bSepDVDL && fplog && do_log)
1596 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1600 /* this factor or 2 correction is necessary
1601 because half of the constraint force is removed
1602 in the vv step, so we have to double it. See
1603 the Redmine issue #1255. It is not yet clear
1604 if the factor of 2 is exact, or just a very
1605 good approximation, and this will be
1606 investigated. The next step is to see if this
1607 can be done adding a dhdl contribution from the
1608 rattle step, but this is somewhat more
1609 complicated with the current code. Will be
1610 investigated, hopefully for 4.6.3. However,
1611 this current solution is much better than
1612 having it completely wrong.
1614 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1618 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1623 /* Need to unshift here */
1624 unshift_self(graph, state->box, state->x);
1629 wallcycle_start(wcycle, ewcVSITECONSTR);
1632 shift_self(graph, state->box, state->x);
1634 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1635 top->idef.iparams, top->idef.il,
1636 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1640 unshift_self(graph, state->box, state->x);
1642 wallcycle_stop(wcycle, ewcVSITECONSTR);
1645 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1646 /* With Leap-Frog we can skip compute_globals at
1647 * non-communication steps, but we need to calculate
1648 * the kinetic energy one step before communication.
1650 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1652 if (ir->nstlist == -1 && bFirstIterate)
1654 gs.sig[eglsNABNSB] = nlh.nabnsb;
1656 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1657 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1659 bFirstIterate ? &gs : NULL,
1660 (step_rel % gs.nstms == 0) &&
1661 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1663 top_global, &bSumEkinhOld,
1665 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1666 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1667 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1668 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1669 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1670 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1673 if (ir->nstlist == -1 && bFirstIterate)
1675 nlh.nabnsb = gs.set[eglsNABNSB];
1676 gs.set[eglsNABNSB] = 0;
1679 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1680 /* ############# END CALC EKIN AND PRESSURE ################# */
1682 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1683 the virial that should probably be addressed eventually. state->veta has better properies,
1684 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1685 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1687 if (iterate.bIterationActive &&
1688 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1689 trace(shake_vir), &tracevir))
1693 bFirstIterate = FALSE;
1696 if (!bVV || bRerunMD)
1698 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1699 sum_dhdl(enerd, state->lambda, ir->fepvals);
1701 update_box(fplog, step, ir, mdatoms, state, f,
1702 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1704 /* ################# END UPDATE STEP 2 ################# */
1705 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1707 /* The coordinates (x) were unshifted in update */
1710 /* We will not sum ekinh_old,
1711 * so signal that we still have to do it.
1713 bSumEkinhOld = TRUE;
1716 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1718 /* use the directly determined last velocity, not actually the averaged half steps */
1719 if (bTrotter && ir->eI == eiVV)
1721 enerd->term[F_EKIN] = last_ekin;
1723 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1727 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1731 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1733 /* ######### END PREPARING EDR OUTPUT ########### */
1738 gmx_bool do_dr, do_or;
1740 if (fplog && do_log && bDoExpanded)
1742 /* only needed if doing expanded ensemble */
1743 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1744 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1746 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1750 upd_mdebin(mdebin, bDoDHDL, TRUE,
1751 t, mdatoms->tmass, enerd, state,
1752 ir->fepvals, ir->expandedvals, lastbox,
1753 shake_vir, force_vir, total_vir, pres,
1754 ekind, mu_tot, constr);
1758 upd_mdebin_step(mdebin);
1761 do_dr = do_per_step(step, ir->nstdisreout);
1762 do_or = do_per_step(step, ir->nstorireout);
1764 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
1766 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1768 if (ir->ePull != epullNO)
1770 pull_print_output(ir->pull, step, t);
1773 if (do_per_step(step, ir->nstlog))
1775 if (fflush(fplog) != 0)
1777 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1783 /* Have to do this part _after_ outputting the logfile and the edr file */
1784 /* Gets written into the state at the beginning of next loop*/
1785 state->fep_state = lamnew;
1787 /* Print the remaining wall clock time for the run */
1788 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1792 fprintf(stderr, "\n");
1794 print_time(stderr, walltime_accounting, step, ir, cr);
1797 /* Replica exchange */
1799 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1800 do_per_step(step, repl_ex_nst))
1802 bExchanged = replica_exchange(fplog, cr, repl_ex,
1803 state_global, enerd,
1806 if (bExchanged && DOMAINDECOMP(cr))
1808 dd_partition_system(fplog, step, cr, TRUE, 1,
1809 state_global, top_global, ir,
1810 state, &f, mdatoms, top, fr,
1811 vsite, shellfc, constr,
1812 nrnb, wcycle, FALSE);
1818 bStartingFromCpt = FALSE;
1820 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1821 /* With all integrators, except VV, we need to retain the pressure
1822 * at the current step for coupling at the next step.
1824 if ((state->flags & (1<<estPRES_PREV)) &&
1826 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1828 /* Store the pressure in t_state for pressure coupling
1829 * at the next MD step.
1831 copy_mat(pres, state->pres_prev);
1834 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1836 if ( (membed != NULL) && (!bLastStep) )
1838 rescale_membed(step_rel, membed, state_global->x);
1845 /* read next frame from input trajectory */
1846 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1851 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1855 if (!bRerunMD || !rerun_fr.bStep)
1857 /* increase the MD step number */
1862 cycles = wallcycle_stop(wcycle, ewcSTEP);
1863 if (DOMAINDECOMP(cr) && wcycle)
1865 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1868 if (bPMETuneRunning || bPMETuneTry)
1870 /* PME grid + cut-off optimization with GPUs or PME nodes */
1872 /* Count the total cycles over the last steps */
1873 cycles_pmes += cycles;
1875 /* We can only switch cut-off at NS steps */
1876 if (step % ir->nstlist == 0)
1878 /* PME grid + cut-off optimization with GPUs or PME nodes */
1881 if (DDMASTER(cr->dd))
1883 /* PME node load is too high, start tuning */
1884 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1886 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1888 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1890 bPMETuneTry = FALSE;
1893 if (bPMETuneRunning)
1895 /* init_step might not be a multiple of nstlist,
1896 * but the first cycle is always skipped anyhow.
1899 pme_load_balance(pme_loadbal, cr,
1900 (bVerbose && MASTER(cr)) ? stderr : NULL,
1902 ir, state, cycles_pmes,
1903 fr->ic, fr->nbv, &fr->pmedata,
1906 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1907 fr->ewaldcoeff = fr->ic->ewaldcoeff;
1908 fr->rlist = fr->ic->rlist;
1909 fr->rlistlong = fr->ic->rlistlong;
1910 fr->rcoulomb = fr->ic->rcoulomb;
1911 fr->rvdw = fr->ic->rvdw;
1917 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1918 gs.set[eglsRESETCOUNTERS] != 0)
1920 /* Reset all the counters related to performance over the run */
1921 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1922 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1923 wcycle_set_reset_counters(wcycle, -1);
1924 if (!(cr->duty & DUTY_PME))
1926 /* Tell our PME node to reset its counters */
1927 gmx_pme_send_resetcounters(cr, step);
1929 /* Correct max_hours for the elapsed time */
1930 max_hours -= elapsed_time/(60.0*60.0);
1931 bResetCountersHalfMaxH = FALSE;
1932 gs.set[eglsRESETCOUNTERS] = 0;
1936 /* End of main MD loop */
1939 /* Stop measuring walltime */
1940 walltime_accounting_end(walltime_accounting);
1942 if (bRerunMD && MASTER(cr))
1947 if (!(cr->duty & DUTY_PME))
1949 /* Tell the PME only node to finish */
1950 gmx_pme_send_finish(cr);
1955 if (ir->nstcalcenergy > 0 && !bRerunMD)
1957 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
1958 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1966 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1968 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)));
1969 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1972 if (pme_loadbal != NULL)
1974 pme_loadbal_done(pme_loadbal, cr, fplog,
1975 fr->nbv != NULL && fr->nbv->bUseGPU);
1978 if (shellfc && fplog)
1980 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1981 (nconverged*100.0)/step_rel);
1982 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1986 if (repl_ex_nst > 0 && MASTER(cr))
1988 print_replica_exchange_statistics(fplog, repl_ex);
1991 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);