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,2014, 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"
68 #include "domdec_network.h"
73 #include "compute_io.h"
75 #include "checkpoint.h"
76 #include "mtop_util.h"
77 #include "sighandler.h"
80 #include "pme_loadbal.h"
83 #include "types/nlistheuristics.h"
84 #include "types/iteratedconstraints.h"
85 #include "nbnxn_cuda_data_mgmt.h"
87 #include "gromacs/utility/gmxmpi.h"
88 #include "gromacs/fileio/confio.h"
89 #include "gromacs/fileio/trajectory_writing.h"
90 #include "gromacs/fileio/trnio.h"
91 #include "gromacs/fileio/trxio.h"
92 #include "gromacs/fileio/xtcio.h"
93 #include "gromacs/timing/wallcycle.h"
94 #include "gromacs/timing/walltime_accounting.h"
95 #include "gromacs/pulling/pull.h"
96 #include "gromacs/swap/swapcoords.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)
152 gmx_mdoutf_t outf = NULL;
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, bNeedRepartition;
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");
289 ir->nstxout_compressed = 0;
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 walltime_accounting_start(walltime_accounting);
672 wallcycle_start(wcycle, ewcRUN);
673 print_start(fplog, cr, walltime_accounting, "mdrun");
675 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
677 chkpt_ret = fcCheckPointParallel( cr->nodeid,
681 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
686 /***********************************************************
690 ************************************************************/
692 /* if rerunMD then read coordinates and velocities from input trajectory */
695 if (getenv("GMX_FORCE_UPDATE"))
703 bNotLastFrame = read_first_frame(oenv, &status,
704 opt2fn("-rerun", nfile, fnm),
705 &rerun_fr, TRX_NEED_X | TRX_READ_V);
706 if (rerun_fr.natoms != top_global->natoms)
709 "Number of atoms in trajectory (%d) does not match the "
710 "run input file (%d)\n",
711 rerun_fr.natoms, top_global->natoms);
713 if (ir->ePBC != epbcNONE)
717 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);
719 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
721 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
728 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
731 if (ir->ePBC != epbcNONE)
733 /* Set the shift vectors.
734 * Necessary here when have a static box different from the tpr box.
736 calc_shifts(rerun_fr.box, fr->shift_vec);
740 /* loop over MD steps or if rerunMD to end of input trajectory */
742 /* Skip the first Nose-Hoover integration when we get the state from tpx */
743 bStateFromTPX = !bStateFromCP;
744 bInitStep = bFirstStep && (bStateFromTPX || bVV);
745 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
747 bSumEkinhOld = FALSE;
749 bNeedRepartition = FALSE;
751 init_global_signals(&gs, cr, ir, repl_ex_nst);
753 step = ir->init_step;
756 if (ir->nstlist == -1)
758 init_nlistheuristics(&nlh, bGStatEveryStep, step);
761 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
763 /* check how many steps are left in other sims */
764 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
768 /* and stop now if we should */
769 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
770 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
771 while (!bLastStep || (bRerunMD && bNotLastFrame))
774 wallcycle_start(wcycle, ewcSTEP);
780 step = rerun_fr.step;
781 step_rel = step - ir->init_step;
794 bLastStep = (step_rel == ir->nsteps);
795 t = t0 + step*ir->delta_t;
798 if (ir->efep != efepNO || ir->bSimTemp)
800 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
801 requiring different logic. */
803 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
804 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
805 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
806 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
807 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
812 update_annealing_target_temp(&(ir->opts), t);
817 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
819 for (i = 0; i < state_global->natoms; i++)
821 copy_rvec(rerun_fr.x[i], state_global->x[i]);
825 for (i = 0; i < state_global->natoms; i++)
827 copy_rvec(rerun_fr.v[i], state_global->v[i]);
832 for (i = 0; i < state_global->natoms; i++)
834 clear_rvec(state_global->v[i]);
838 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
839 " Ekin, temperature and pressure are incorrect,\n"
840 " the virial will be incorrect when constraints are present.\n"
842 bRerunWarnNoV = FALSE;
846 copy_mat(rerun_fr.box, state_global->box);
847 copy_mat(state_global->box, state->box);
849 if (vsite && (Flags & MD_RERUN_VSITE))
851 if (DOMAINDECOMP(cr))
853 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
857 /* Following is necessary because the graph may get out of sync
858 * with the coordinates if we only have every N'th coordinate set
860 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
861 shift_self(graph, state->box, state->x);
863 construct_vsites(vsite, state->x, ir->delta_t, state->v,
864 top->idef.iparams, top->idef.il,
865 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
868 unshift_self(graph, state->box, state->x);
873 /* Stop Center of Mass motion */
874 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
878 /* for rerun MD always do Neighbour Searching */
879 bNS = (bFirstStep || ir->nstlist != 0);
884 /* Determine whether or not to do Neighbour Searching and LR */
885 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
887 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
888 (ir->nstlist == -1 && nlh.nabnsb > 0));
890 if (bNS && ir->nstlist == -1)
892 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
896 /* check whether we should stop because another simulation has
900 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
901 (multisim_nsteps != ir->nsteps) )
908 "Stopping simulation %d because another one has finished\n",
912 gs.sig[eglsCHKPT] = 1;
917 /* < 0 means stop at next step, > 0 means stop at next NS step */
918 if ( (gs.set[eglsSTOPCOND] < 0) ||
919 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
924 /* Determine whether or not to update the Born radii if doing GB */
925 bBornRadii = bFirstStep;
926 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
931 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
932 do_verbose = bVerbose &&
933 (step % stepout == 0 || bFirstStep || bLastStep);
935 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
943 bMasterState = FALSE;
944 /* Correct the new box if it is too skewed */
945 if (DYNAMIC_BOX(*ir))
947 if (correct_box(fplog, step, state->box, graph))
952 if (DOMAINDECOMP(cr) && bMasterState)
954 dd_collect_state(cr->dd, state, state_global);
958 if (DOMAINDECOMP(cr))
960 /* Repartition the domain decomposition */
961 wallcycle_start(wcycle, ewcDOMDEC);
962 dd_partition_system(fplog, step, cr,
963 bMasterState, nstglobalcomm,
964 state_global, top_global, ir,
965 state, &f, mdatoms, top, fr,
966 vsite, shellfc, constr,
968 do_verbose && !bPMETuneRunning);
969 wallcycle_stop(wcycle, ewcDOMDEC);
970 /* If using an iterative integrator, reallocate space to match the decomposition */
974 if (MASTER(cr) && do_log)
976 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
979 if (ir->efep != efepNO)
981 update_mdatoms(mdatoms, state->lambda[efptMASS]);
984 if ((bRerunMD && rerun_fr.bV) || bExchanged)
987 /* We need the kinetic energy at minus the half step for determining
988 * the full step kinetic energy and possibly for T-coupling.*/
989 /* This may not be quite working correctly yet . . . . */
990 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
991 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
992 constr, NULL, FALSE, state->box,
993 top_global, &bSumEkinhOld,
994 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
996 clear_mat(force_vir);
998 /* We write a checkpoint at this MD step when:
999 * either at an NS step when we signalled through gs,
1000 * or at the last step (but not when we do not want confout),
1001 * but never at the first step or with rerun.
1003 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1004 (bLastStep && (Flags & MD_CONFOUT))) &&
1005 step > ir->init_step && !bRerunMD);
1008 gs.set[eglsCHKPT] = 0;
1011 /* Determine the energy and pressure:
1012 * at nstcalcenergy steps and at energy output steps (set below).
1014 if (EI_VV(ir->eI) && (!bInitStep))
1016 /* for vv, the first half of the integration actually corresponds
1017 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1018 but the virial needs to be calculated on both the current step and the 'next' step. Future
1019 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1021 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1022 bCalcVir = bCalcEner ||
1023 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1027 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1028 bCalcVir = bCalcEner ||
1029 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1032 /* Do we need global communication ? */
1033 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1034 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1035 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1037 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1039 if (do_ene || do_log)
1046 /* these CGLO_ options remain the same throughout the iteration */
1047 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1048 (bGStat ? CGLO_GSTAT : 0)
1051 force_flags = (GMX_FORCE_STATECHANGED |
1052 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1053 GMX_FORCE_ALLFORCES |
1055 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1056 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1057 (bDoFEP ? GMX_FORCE_DHDL : 0)
1062 if (do_per_step(step, ir->nstcalclr))
1064 force_flags |= GMX_FORCE_DO_LR;
1070 /* Now is the time to relax the shells */
1071 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1072 ir, bNS, force_flags,
1075 state, f, force_vir, mdatoms,
1076 nrnb, wcycle, graph, groups,
1077 shellfc, fr, bBornRadii, t, mu_tot,
1079 mdoutf_get_fp_field(outf));
1089 /* The coordinates (x) are shifted (to get whole molecules)
1091 * This is parallellized as well, and does communication too.
1092 * Check comments in sim_util.c
1094 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1095 state->box, state->x, &state->hist,
1096 f, force_vir, mdatoms, enerd, fcd,
1097 state->lambda, graph,
1098 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1099 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1102 if (bVV && !bStartingFromCpt && !bRerunMD)
1103 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1105 if (ir->eI == eiVV && bInitStep)
1107 /* if using velocity verlet with full time step Ekin,
1108 * take the first half step only to compute the
1109 * virial for the first step. From there,
1110 * revert back to the initial coordinates
1111 * so that the input is actually the initial step.
1113 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1117 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1118 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1121 /* If we are using twin-range interactions where the long-range component
1122 * is only evaluated every nstcalclr>1 steps, we should do a special update
1123 * step to combine the long-range forces on these steps.
1124 * For nstcalclr=1 this is not done, since the forces would have been added
1125 * directly to the short-range forces already.
1127 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1129 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1130 f, bUpdateDoLR, fr->f_twin, fcd,
1131 ekind, M, upd, bInitStep, etrtVELOCITY1,
1132 cr, nrnb, constr, &top->idef);
1134 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1136 gmx_iterate_init(&iterate, TRUE);
1138 /* for iterations, we save these vectors, as we will be self-consistently iterating
1141 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1143 /* save the state */
1144 if (iterate.bIterationActive)
1146 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1149 bFirstIterate = TRUE;
1150 while (bFirstIterate || iterate.bIterationActive)
1152 if (iterate.bIterationActive)
1154 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1155 if (bFirstIterate && bTrotter)
1157 /* The first time through, we need a decent first estimate
1158 of veta(t+dt) to compute the constraints. Do
1159 this by computing the box volume part of the
1160 trotter integration at this time. Nothing else
1161 should be changed by this routine here. If
1162 !(first time), we start with the previous value
1165 veta_save = state->veta;
1166 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1167 vetanew = state->veta;
1168 state->veta = veta_save;
1173 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1175 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1176 state, fr->bMolPBC, graph, f,
1177 &top->idef, shake_vir,
1178 cr, nrnb, wcycle, upd, constr,
1179 TRUE, bCalcVir, vetanew);
1183 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1189 /* Need to unshift here if a do_force has been
1190 called in the previous step */
1191 unshift_self(graph, state->box, state->x);
1194 /* if VV, compute the pressure and constraints */
1195 /* For VV2, we strictly only need this if using pressure
1196 * control, but we really would like to have accurate pressures
1198 * Think about ways around this in the future?
1199 * For now, keep this choice in comments.
1201 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1202 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1204 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1205 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1207 bSumEkinhOld = TRUE;
1209 /* for vv, the first half of the integration actually corresponds to the previous step.
1210 So we need information from the last step in the first half of the integration */
1211 if (bGStat || do_per_step(step-1, nstglobalcomm))
1213 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1214 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1215 constr, NULL, FALSE, state->box,
1216 top_global, &bSumEkinhOld,
1219 | (bTemp ? CGLO_TEMPERATURE : 0)
1220 | (bPres ? CGLO_PRESSURE : 0)
1221 | (bPres ? CGLO_CONSTRAINT : 0)
1222 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1223 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1226 /* explanation of above:
1227 a) We compute Ekin at the full time step
1228 if 1) we are using the AveVel Ekin, and it's not the
1229 initial step, or 2) if we are using AveEkin, but need the full
1230 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1231 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1232 EkinAveVel because it's needed for the pressure */
1234 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1239 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1240 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1247 /* We need the kinetic energy at minus the half step for determining
1248 * the full step kinetic energy and possibly for T-coupling.*/
1249 /* This may not be quite working correctly yet . . . . */
1250 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1251 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1252 constr, NULL, FALSE, state->box,
1253 top_global, &bSumEkinhOld,
1254 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1259 if (iterate.bIterationActive &&
1260 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1261 state->veta, &vetanew))
1265 bFirstIterate = FALSE;
1268 if (bTrotter && !bInitStep)
1270 copy_mat(shake_vir, state->svir_prev);
1271 copy_mat(force_vir, state->fvir_prev);
1272 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1274 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1275 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1276 enerd->term[F_EKIN] = trace(ekind->ekin);
1279 /* if it's the initial step, we performed this first step just to get the constraint virial */
1280 if (bInitStep && ir->eI == eiVV)
1282 copy_rvecn(cbuf, state->v, 0, state->natoms);
1286 /* MRS -- now done iterating -- compute the conserved quantity */
1289 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1292 last_ekin = enerd->term[F_EKIN];
1294 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1296 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1298 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1301 sum_dhdl(enerd, state->lambda, ir->fepvals);
1305 /* ######## END FIRST UPDATE STEP ############## */
1306 /* ######## If doing VV, we now have v(dt) ###### */
1309 /* perform extended ensemble sampling in lambda - we don't
1310 actually move to the new state before outputting
1311 statistics, but if performing simulated tempering, we
1312 do update the velocities and the tau_t. */
1314 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, mcrng, state->v, mdatoms);
1315 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1316 copy_df_history(&state_global->dfhist, &state->dfhist);
1319 /* Now we have the energies and forces corresponding to the
1320 * coordinates at time t. We must output all of this before
1323 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1324 ir, state, state_global, top_global, fr, upd,
1325 outf, mdebin, ekind, f, f_global,
1326 wcycle, mcrng, &nchkpt,
1327 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1330 /* kludge -- virial is lost with restart for NPT control. Must restart */
1331 if (bStartingFromCpt && bVV)
1333 copy_mat(state->svir_prev, shake_vir);
1334 copy_mat(state->fvir_prev, force_vir);
1337 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1339 /* Check whether everything is still allright */
1340 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1341 #ifdef GMX_THREAD_MPI
1346 /* this is just make gs.sig compatible with the hack
1347 of sending signals around by MPI_Reduce with together with
1349 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1351 gs.sig[eglsSTOPCOND] = 1;
1353 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1355 gs.sig[eglsSTOPCOND] = -1;
1357 /* < 0 means stop at next step, > 0 means stop at next NS step */
1361 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1362 gmx_get_signal_name(),
1363 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1367 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1368 gmx_get_signal_name(),
1369 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1371 handled_stop_condition = (int)gmx_get_stop_condition();
1373 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1374 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1375 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1377 /* Signal to terminate the run */
1378 gs.sig[eglsSTOPCOND] = 1;
1381 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1383 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1386 if (bResetCountersHalfMaxH && MASTER(cr) &&
1387 elapsed_time > max_hours*60.0*60.0*0.495)
1389 gs.sig[eglsRESETCOUNTERS] = 1;
1392 if (ir->nstlist == -1 && !bRerunMD)
1394 /* When bGStatEveryStep=FALSE, global_stat is only called
1395 * when we check the atom displacements, not at NS steps.
1396 * This means that also the bonded interaction count check is not
1397 * performed immediately after NS. Therefore a few MD steps could
1398 * be performed with missing interactions.
1399 * But wrong energies are never written to file,
1400 * since energies are only written after global_stat
1403 if (step >= nlh.step_nscheck)
1405 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1406 nlh.scale_tot, state->x);
1410 /* This is not necessarily true,
1411 * but step_nscheck is determined quite conservatively.
1417 /* In parallel we only have to check for checkpointing in steps
1418 * where we do global communication,
1419 * otherwise the other nodes don't know.
1421 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1424 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1425 gs.set[eglsCHKPT] == 0)
1427 gs.sig[eglsCHKPT] = 1;
1430 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1435 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1437 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1439 gmx_bool bIfRandomize;
1440 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1441 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1442 if (constr && bIfRandomize)
1444 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1445 state, fr->bMolPBC, graph, f,
1446 &top->idef, tmp_vir,
1447 cr, nrnb, wcycle, upd, constr,
1448 TRUE, bCalcVir, vetanew);
1453 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1455 gmx_iterate_init(&iterate, TRUE);
1456 /* for iterations, we save these vectors, as we will be redoing the calculations */
1457 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1460 bFirstIterate = TRUE;
1461 while (bFirstIterate || iterate.bIterationActive)
1463 /* We now restore these vectors to redo the calculation with improved extended variables */
1464 if (iterate.bIterationActive)
1466 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1469 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1470 so scroll down for that logic */
1472 /* ######### START SECOND UPDATE STEP ################# */
1473 /* Box is changed in update() when we do pressure coupling,
1474 * but we should still use the old box for energy corrections and when
1475 * writing it to the energy file, so it matches the trajectory files for
1476 * the same timestep above. Make a copy in a separate array.
1478 copy_mat(state->box, lastbox);
1483 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1485 wallcycle_start(wcycle, ewcUPDATE);
1486 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1489 if (iterate.bIterationActive)
1497 /* we use a new value of scalevir to converge the iterations faster */
1498 scalevir = tracevir/trace(shake_vir);
1500 msmul(shake_vir, scalevir, shake_vir);
1501 m_add(force_vir, shake_vir, total_vir);
1502 clear_mat(shake_vir);
1504 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1505 /* We can only do Berendsen coupling after we have summed
1506 * the kinetic energy or virial. Since the happens
1507 * in global_state after update, we should only do it at
1508 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1513 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1514 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1519 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1521 /* velocity half-step update */
1522 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1523 bUpdateDoLR, fr->f_twin, fcd,
1524 ekind, M, upd, FALSE, etrtVELOCITY2,
1525 cr, nrnb, constr, &top->idef);
1528 /* Above, initialize just copies ekinh into ekin,
1529 * it doesn't copy position (for VV),
1530 * and entire integrator for MD.
1533 if (ir->eI == eiVVAK)
1535 copy_rvecn(state->x, cbuf, 0, state->natoms);
1537 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1539 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1540 bUpdateDoLR, fr->f_twin, fcd,
1541 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1542 wallcycle_stop(wcycle, ewcUPDATE);
1544 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1545 fr->bMolPBC, graph, f,
1546 &top->idef, shake_vir,
1547 cr, nrnb, wcycle, upd, constr,
1548 FALSE, bCalcVir, state->veta);
1550 if (ir->eI == eiVVAK)
1552 /* erase F_EKIN and F_TEMP here? */
1553 /* just compute the kinetic energy at the half step to perform a trotter step */
1554 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1555 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1556 constr, NULL, FALSE, lastbox,
1557 top_global, &bSumEkinhOld,
1558 cglo_flags | CGLO_TEMPERATURE
1560 wallcycle_start(wcycle, ewcUPDATE);
1561 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1562 /* now we know the scaling, we can compute the positions again again */
1563 copy_rvecn(cbuf, state->x, 0, state->natoms);
1565 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1567 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1568 bUpdateDoLR, fr->f_twin, fcd,
1569 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1570 wallcycle_stop(wcycle, ewcUPDATE);
1572 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1573 /* are the small terms in the shake_vir here due
1574 * to numerical errors, or are they important
1575 * physically? I'm thinking they are just errors, but not completely sure.
1576 * For now, will call without actually constraining, constr=NULL*/
1577 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1578 state, fr->bMolPBC, graph, f,
1579 &top->idef, tmp_vir,
1580 cr, nrnb, wcycle, upd, NULL,
1586 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1589 if (fr->bSepDVDL && fplog && do_log)
1591 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1595 /* this factor or 2 correction is necessary
1596 because half of the constraint force is removed
1597 in the vv step, so we have to double it. See
1598 the Redmine issue #1255. It is not yet clear
1599 if the factor of 2 is exact, or just a very
1600 good approximation, and this will be
1601 investigated. The next step is to see if this
1602 can be done adding a dhdl contribution from the
1603 rattle step, but this is somewhat more
1604 complicated with the current code. Will be
1605 investigated, hopefully for 4.6.3. However,
1606 this current solution is much better than
1607 having it completely wrong.
1609 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1613 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1618 /* Need to unshift here */
1619 unshift_self(graph, state->box, state->x);
1624 wallcycle_start(wcycle, ewcVSITECONSTR);
1627 shift_self(graph, state->box, state->x);
1629 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1630 top->idef.iparams, top->idef.il,
1631 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1635 unshift_self(graph, state->box, state->x);
1637 wallcycle_stop(wcycle, ewcVSITECONSTR);
1640 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1641 /* With Leap-Frog we can skip compute_globals at
1642 * non-communication steps, but we need to calculate
1643 * the kinetic energy one step before communication.
1645 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1647 if (ir->nstlist == -1 && bFirstIterate)
1649 gs.sig[eglsNABNSB] = nlh.nabnsb;
1651 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1652 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1654 bFirstIterate ? &gs : NULL,
1655 (step_rel % gs.nstms == 0) &&
1656 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1658 top_global, &bSumEkinhOld,
1660 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1661 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1662 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1663 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1664 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1665 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1668 if (ir->nstlist == -1 && bFirstIterate)
1670 nlh.nabnsb = gs.set[eglsNABNSB];
1671 gs.set[eglsNABNSB] = 0;
1674 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1675 /* ############# END CALC EKIN AND PRESSURE ################# */
1677 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1678 the virial that should probably be addressed eventually. state->veta has better properies,
1679 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1680 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1682 if (iterate.bIterationActive &&
1683 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1684 trace(shake_vir), &tracevir))
1688 bFirstIterate = FALSE;
1691 if (!bVV || bRerunMD)
1693 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1694 sum_dhdl(enerd, state->lambda, ir->fepvals);
1696 update_box(fplog, step, ir, mdatoms, state, f,
1697 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1699 /* ################# END UPDATE STEP 2 ################# */
1700 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1702 /* The coordinates (x) were unshifted in update */
1705 /* We will not sum ekinh_old,
1706 * so signal that we still have to do it.
1708 bSumEkinhOld = TRUE;
1711 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1713 /* use the directly determined last velocity, not actually the averaged half steps */
1714 if (bTrotter && ir->eI == eiVV)
1716 enerd->term[F_EKIN] = last_ekin;
1718 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1722 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1726 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1728 /* ######### END PREPARING EDR OUTPUT ########### */
1733 gmx_bool do_dr, do_or;
1735 if (fplog && do_log && bDoExpanded)
1737 /* only needed if doing expanded ensemble */
1738 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1739 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1741 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1745 upd_mdebin(mdebin, bDoDHDL, TRUE,
1746 t, mdatoms->tmass, enerd, state,
1747 ir->fepvals, ir->expandedvals, lastbox,
1748 shake_vir, force_vir, total_vir, pres,
1749 ekind, mu_tot, constr);
1753 upd_mdebin_step(mdebin);
1756 do_dr = do_per_step(step, ir->nstdisreout);
1757 do_or = do_per_step(step, ir->nstorireout);
1759 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1761 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1763 if (ir->ePull != epullNO)
1765 pull_print_output(ir->pull, step, t);
1768 if (do_per_step(step, ir->nstlog))
1770 if (fflush(fplog) != 0)
1772 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1778 /* Have to do this part _after_ outputting the logfile and the edr file */
1779 /* Gets written into the state at the beginning of next loop*/
1780 state->fep_state = lamnew;
1782 /* Print the remaining wall clock time for the run */
1783 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1787 fprintf(stderr, "\n");
1789 print_time(stderr, walltime_accounting, step, ir, cr);
1792 /* Ion/water position swapping.
1793 * Not done in last step since trajectory writing happens before this call
1794 * in the MD loop and exchanges would be lost anyway. */
1795 bNeedRepartition = FALSE;
1796 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1797 do_per_step(step, ir->swap->nstswap))
1799 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1800 bRerunMD ? rerun_fr.x : state->x,
1801 bRerunMD ? rerun_fr.box : state->box,
1802 top_global, MASTER(cr) && bVerbose, bRerunMD);
1804 if (bNeedRepartition && DOMAINDECOMP(cr))
1806 dd_collect_state(cr->dd, state, state_global);
1810 /* Replica exchange */
1812 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1813 do_per_step(step, repl_ex_nst))
1815 bExchanged = replica_exchange(fplog, cr, repl_ex,
1816 state_global, enerd,
1820 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1822 dd_partition_system(fplog, step, cr, TRUE, 1,
1823 state_global, top_global, ir,
1824 state, &f, mdatoms, top, fr,
1825 vsite, shellfc, constr,
1826 nrnb, wcycle, FALSE);
1831 bStartingFromCpt = FALSE;
1833 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1834 /* With all integrators, except VV, we need to retain the pressure
1835 * at the current step for coupling at the next step.
1837 if ((state->flags & (1<<estPRES_PREV)) &&
1839 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1841 /* Store the pressure in t_state for pressure coupling
1842 * at the next MD step.
1844 copy_mat(pres, state->pres_prev);
1847 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1849 if ( (membed != NULL) && (!bLastStep) )
1851 rescale_membed(step_rel, membed, state_global->x);
1858 /* read next frame from input trajectory */
1859 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1864 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1868 if (!bRerunMD || !rerun_fr.bStep)
1870 /* increase the MD step number */
1875 cycles = wallcycle_stop(wcycle, ewcSTEP);
1876 if (DOMAINDECOMP(cr) && wcycle)
1878 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1881 if (bPMETuneRunning || bPMETuneTry)
1883 /* PME grid + cut-off optimization with GPUs or PME nodes */
1885 /* Count the total cycles over the last steps */
1886 cycles_pmes += cycles;
1888 /* We can only switch cut-off at NS steps */
1889 if (step % ir->nstlist == 0)
1891 /* PME grid + cut-off optimization with GPUs or PME nodes */
1894 if (DDMASTER(cr->dd))
1896 /* PME node load is too high, start tuning */
1897 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1899 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1901 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1903 bPMETuneTry = FALSE;
1906 if (bPMETuneRunning)
1908 /* init_step might not be a multiple of nstlist,
1909 * but the first cycle is always skipped anyhow.
1912 pme_load_balance(pme_loadbal, cr,
1913 (bVerbose && MASTER(cr)) ? stderr : NULL,
1915 ir, state, cycles_pmes,
1916 fr->ic, fr->nbv, &fr->pmedata,
1919 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1920 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1921 fr->rlist = fr->ic->rlist;
1922 fr->rlistlong = fr->ic->rlistlong;
1923 fr->rcoulomb = fr->ic->rcoulomb;
1924 fr->rvdw = fr->ic->rvdw;
1930 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1931 gs.set[eglsRESETCOUNTERS] != 0)
1933 /* Reset all the counters related to performance over the run */
1934 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1935 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1936 wcycle_set_reset_counters(wcycle, -1);
1937 if (!(cr->duty & DUTY_PME))
1939 /* Tell our PME node to reset its counters */
1940 gmx_pme_send_resetcounters(cr, step);
1942 /* Correct max_hours for the elapsed time */
1943 max_hours -= elapsed_time/(60.0*60.0);
1944 bResetCountersHalfMaxH = FALSE;
1945 gs.set[eglsRESETCOUNTERS] = 0;
1949 /* End of main MD loop */
1952 /* Stop measuring walltime */
1953 walltime_accounting_end(walltime_accounting);
1955 if (bRerunMD && MASTER(cr))
1960 if (!(cr->duty & DUTY_PME))
1962 /* Tell the PME only node to finish */
1963 gmx_pme_send_finish(cr);
1968 if (ir->nstcalcenergy > 0 && !bRerunMD)
1970 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1971 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1978 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1980 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)));
1981 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1984 if (pme_loadbal != NULL)
1986 pme_loadbal_done(pme_loadbal, cr, fplog,
1987 fr->nbv != NULL && fr->nbv->bUseGPU);
1990 if (shellfc && fplog)
1992 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1993 (nconverged*100.0)/step_rel);
1994 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1998 if (repl_ex_nst > 0 && MASTER(cr))
2000 print_replica_exchange_statistics(fplog, repl_ex);
2003 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);