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,2015, 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.
42 #include "gromacs/utility/smalloc.h"
54 #include "md_support.h"
55 #include "md_logging.h"
69 #include "domdec_network.h"
70 #include "gromacs/gmxlib/topsort.h"
74 #include "gromacs/gmxpreprocess/compute_io.h"
75 #include "checkpoint.h"
76 #include "mtop_util.h"
77 #include "sighandler.h"
79 #include "gromacs/utility/cstringutil.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"
97 #include "gromacs/imd/imd.h"
100 #include "corewrap.h"
103 static void reset_all_counters(FILE *fplog, t_commrec *cr,
105 gmx_int64_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", gmx_gettime());
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,
152 gmx_walltime_accounting_t walltime_accounting)
154 gmx_mdoutf_t outf = NULL;
155 gmx_int64_t step, step_rel;
157 double t, t0, lam0[efptNR];
158 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
159 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
160 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
161 bBornRadii, bStartingFromCpt;
162 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
163 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
164 bForceUpdate = FALSE, bCPT;
165 gmx_bool bMasterState;
166 int force_flags, cglo_flags;
167 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
172 t_state *bufstate = NULL;
173 matrix *scale_tot, pcoupl_mu, M, ebox;
176 gmx_repl_ex_t repl_ex = NULL;
179 t_mdebin *mdebin = NULL;
180 t_state *state = NULL;
181 rvec *f_global = NULL;
182 gmx_enerdata_t *enerd;
184 gmx_global_stat_t gstat;
185 gmx_update_t upd = NULL;
186 t_graph *graph = 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, bDoReplEx, bExchanged, bNeedRepartition;
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_int64_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 gmx_bool bIMDstep = FALSE;
229 /* Temporary addition for FAHCORE checkpointing */
233 /* Check for special mdrun options */
234 bRerunMD = (Flags & MD_RERUN);
235 bAppend = (Flags & MD_APPENDFILES);
236 if (Flags & MD_RESETCOUNTERSHALFWAY)
240 /* Signal to reset the counters half the simulation steps. */
241 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
243 /* Signal to reset the counters halfway the simulation time. */
244 bResetCountersHalfMaxH = (max_hours > 0);
247 /* md-vv uses averaged full step velocities for T-control
248 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
249 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
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, wcycle);
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);
333 if (shellfc && ir->nstcalcenergy != 1)
335 gmx_fatal(FARGS, "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is not supported in combinations with shell particles.\nPlease make a new tpr file.", ir->nstcalcenergy);
337 if (shellfc && DOMAINDECOMP(cr))
339 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
341 if (shellfc && ir->eI == eiNM)
343 /* Currently shells don't work with Normal Modes */
344 gmx_fatal(FARGS, "Normal Mode analysis is not supported with shells.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
347 if (vsite && ir->eI == eiNM)
349 /* Currently virtual sites don't work with Normal Modes */
350 gmx_fatal(FARGS, "Normal Mode analysis is not supported with virtual sites.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
355 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
356 set_deform_reference_box(upd,
357 deform_init_init_step_tpx,
358 deform_init_box_tpx);
359 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
363 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
364 if ((io > 2000) && MASTER(cr))
367 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
372 if (DOMAINDECOMP(cr))
374 top = dd_init_local_top(top_global);
377 dd_init_local_state(cr->dd, state_global, state);
379 if (DDMASTER(cr->dd) && ir->nstfout)
381 snew(f_global, state_global->natoms);
386 top = gmx_mtop_generate_local_top(top_global, ir);
388 forcerec_set_excl_load(fr, top);
390 state = serial_init_local_state(state_global);
393 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
397 set_vsite_top(vsite, top, mdatoms, cr);
400 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
402 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
407 make_local_shells(cr, mdatoms, shellfc);
410 setup_bonded_threading(fr, &top->idef);
413 /* Set up interactive MD (IMD) */
414 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
415 nfile, fnm, oenv, imdport, Flags);
417 if (DOMAINDECOMP(cr))
419 /* Distribute the charge groups over the nodes from the master node */
420 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
421 state_global, top_global, ir,
422 state, &f, mdatoms, top, fr,
423 vsite, shellfc, constr,
424 nrnb, wcycle, FALSE);
428 update_mdatoms(mdatoms, state->lambda[efptMASS]);
430 if (opt2bSet("-cpi", nfile, fnm))
432 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
436 bStateFromCP = FALSE;
441 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
448 /* Update mdebin with energy history if appending to output files */
449 if (Flags & MD_APPENDFILES)
451 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
455 /* We might have read an energy history from checkpoint,
456 * free the allocated memory and reset the counts.
458 done_energyhistory(&state_global->enerhist);
459 init_energyhistory(&state_global->enerhist);
462 /* Set the initial energy history in state by updating once */
463 update_energyhistory(&state_global->enerhist, mdebin);
466 /* Initialize constraints */
467 if (constr && !DOMAINDECOMP(cr))
469 set_constraints(constr, top, ir, mdatoms, cr);
472 if (repl_ex_nst > 0 && MASTER(cr))
474 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
475 repl_ex_nst, repl_ex_nex, repl_ex_seed);
478 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
479 * PME tuning is not supported with PME only for LJ and not for Coulomb.
481 if ((Flags & MD_TUNEPME) &&
482 EEL_PME(fr->eeltype) &&
483 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
486 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
488 if (cr->duty & DUTY_PME)
490 /* Start tuning right away, as we can't measure the load */
491 bPMETuneRunning = TRUE;
495 /* Separate PME nodes, we can measure the PP/PME load balance */
500 if (!ir->bContinuation && !bRerunMD)
502 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
504 /* Set the velocities of frozen particles to zero */
505 for (i = 0; i < mdatoms->homenr; i++)
507 for (m = 0; m < DIM; m++)
509 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
519 /* Constrain the initial coordinates and velocities */
520 do_constrain_first(fplog, constr, ir, mdatoms, state,
525 /* Construct the virtual sites for the initial configuration */
526 construct_vsites(vsite, state->x, ir->delta_t, NULL,
527 top->idef.iparams, top->idef.il,
528 fr->ePBC, fr->bMolPBC, cr, state->box);
534 /* set free energy calculation frequency as the minimum
535 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
536 nstfep = ir->fepvals->nstdhdl;
539 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
543 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
546 /* I'm assuming we need global communication the first time! MRS */
547 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
548 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
549 | (bVV ? CGLO_PRESSURE : 0)
550 | (bVV ? CGLO_CONSTRAINT : 0)
551 | (bRerunMD ? CGLO_RERUNMD : 0)
552 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
554 bSumEkinhOld = FALSE;
555 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
556 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
557 constr, NULL, FALSE, state->box,
558 top_global, &bSumEkinhOld, cglo_flags);
559 if (ir->eI == eiVVAK)
561 /* a second call to get the half step temperature initialized as well */
562 /* we do the same call as above, but turn the pressure off -- internally to
563 compute_globals, this is recognized as a velocity verlet half-step
564 kinetic energy calculation. This minimized excess variables, but
565 perhaps loses some logic?*/
567 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
568 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
569 constr, NULL, FALSE, state->box,
570 top_global, &bSumEkinhOld,
571 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
574 /* Calculate the initial half step temperature, and save the ekinh_old */
575 if (!(Flags & MD_STARTFROMCPT))
577 for (i = 0; (i < ir->opts.ngtc); i++)
579 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
584 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
585 and there is no previous step */
588 /* if using an iterative algorithm, we need to create a working directory for the state. */
591 bufstate = init_bufstate(state);
594 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
595 temperature control */
596 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
600 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
603 "RMS relative constraint deviation after constraining: %.2e\n",
604 constr_rmsd(constr, FALSE));
606 if (EI_STATE_VELOCITY(ir->eI))
608 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
612 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
613 " input trajectory '%s'\n\n",
614 *(top_global->name), opt2fn("-rerun", nfile, fnm));
617 fprintf(stderr, "Calculated time to finish depends on nsteps from "
618 "run input file,\nwhich may not correspond to the time "
619 "needed to process input trajectory.\n\n");
625 fprintf(stderr, "starting mdrun '%s'\n",
626 *(top_global->name));
629 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
633 sprintf(tbuf, "%s", "infinite");
635 if (ir->init_step > 0)
637 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
638 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
639 gmx_step_str(ir->init_step, sbuf2),
640 ir->init_step*ir->delta_t);
644 fprintf(stderr, "%s steps, %s ps.\n",
645 gmx_step_str(ir->nsteps, sbuf), tbuf);
648 fprintf(fplog, "\n");
651 walltime_accounting_start(walltime_accounting);
652 wallcycle_start(wcycle, ewcRUN);
653 print_start(fplog, cr, walltime_accounting, "mdrun");
655 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
657 chkpt_ret = fcCheckPointParallel( cr->nodeid,
661 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
666 /***********************************************************
670 ************************************************************/
672 /* if rerunMD then read coordinates and velocities from input trajectory */
675 if (getenv("GMX_FORCE_UPDATE"))
683 bNotLastFrame = read_first_frame(oenv, &status,
684 opt2fn("-rerun", nfile, fnm),
685 &rerun_fr, TRX_NEED_X | TRX_READ_V);
686 if (rerun_fr.natoms != top_global->natoms)
689 "Number of atoms in trajectory (%d) does not match the "
690 "run input file (%d)\n",
691 rerun_fr.natoms, top_global->natoms);
693 if (ir->ePBC != epbcNONE)
697 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);
699 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
701 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
708 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
711 if (ir->ePBC != epbcNONE)
713 /* Set the shift vectors.
714 * Necessary here when have a static box different from the tpr box.
716 calc_shifts(rerun_fr.box, fr->shift_vec);
720 /* loop over MD steps or if rerunMD to end of input trajectory */
722 /* Skip the first Nose-Hoover integration when we get the state from tpx */
723 bStateFromTPX = !bStateFromCP;
724 bInitStep = bFirstStep && (bStateFromTPX || bVV);
725 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
727 bSumEkinhOld = FALSE;
730 bNeedRepartition = FALSE;
732 init_global_signals(&gs, cr, ir, repl_ex_nst);
734 step = ir->init_step;
737 if (ir->nstlist == -1)
739 init_nlistheuristics(&nlh, bGStatEveryStep, step);
742 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
744 /* check how many steps are left in other sims */
745 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
749 /* and stop now if we should */
750 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
751 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
752 while (!bLastStep || (bRerunMD && bNotLastFrame))
755 wallcycle_start(wcycle, ewcSTEP);
761 step = rerun_fr.step;
762 step_rel = step - ir->init_step;
775 bLastStep = (step_rel == ir->nsteps);
776 t = t0 + step*ir->delta_t;
779 if (ir->efep != efepNO || ir->bSimTemp)
781 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
782 requiring different logic. */
784 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
785 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
786 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
787 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
788 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
791 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
792 do_per_step(step, repl_ex_nst));
796 update_annealing_target_temp(&(ir->opts), t);
801 if (!DOMAINDECOMP(cr) || MASTER(cr))
803 for (i = 0; i < state_global->natoms; i++)
805 copy_rvec(rerun_fr.x[i], state_global->x[i]);
809 for (i = 0; i < state_global->natoms; i++)
811 copy_rvec(rerun_fr.v[i], state_global->v[i]);
816 for (i = 0; i < state_global->natoms; i++)
818 clear_rvec(state_global->v[i]);
822 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
823 " Ekin, temperature and pressure are incorrect,\n"
824 " the virial will be incorrect when constraints are present.\n"
826 bRerunWarnNoV = FALSE;
830 copy_mat(rerun_fr.box, state_global->box);
831 copy_mat(state_global->box, state->box);
833 if (vsite && (Flags & MD_RERUN_VSITE))
835 if (DOMAINDECOMP(cr))
837 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
841 /* Following is necessary because the graph may get out of sync
842 * with the coordinates if we only have every N'th coordinate set
844 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
845 shift_self(graph, state->box, state->x);
847 construct_vsites(vsite, state->x, ir->delta_t, state->v,
848 top->idef.iparams, top->idef.il,
849 fr->ePBC, fr->bMolPBC, cr, state->box);
852 unshift_self(graph, state->box, state->x);
857 /* Stop Center of Mass motion */
858 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
862 /* for rerun MD always do Neighbour Searching */
863 bNS = (bFirstStep || ir->nstlist != 0);
868 /* Determine whether or not to do Neighbour Searching and LR */
869 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
871 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
872 (ir->nstlist == -1 && nlh.nabnsb > 0));
874 if (bNS && ir->nstlist == -1)
876 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
880 /* check whether we should stop because another simulation has
884 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
885 (multisim_nsteps != ir->nsteps) )
892 "Stopping simulation %d because another one has finished\n",
896 gs.sig[eglsCHKPT] = 1;
901 /* < 0 means stop at next step, > 0 means stop at next NS step */
902 if ( (gs.set[eglsSTOPCOND] < 0) ||
903 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
908 /* Determine whether or not to update the Born radii if doing GB */
909 bBornRadii = bFirstStep;
910 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
915 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
916 do_verbose = bVerbose &&
917 (step % stepout == 0 || bFirstStep || bLastStep);
919 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
927 bMasterState = FALSE;
928 /* Correct the new box if it is too skewed */
929 if (DYNAMIC_BOX(*ir))
931 if (correct_box(fplog, step, state->box, graph))
936 if (DOMAINDECOMP(cr) && bMasterState)
938 dd_collect_state(cr->dd, state, state_global);
942 if (DOMAINDECOMP(cr))
944 /* Repartition the domain decomposition */
945 wallcycle_start(wcycle, ewcDOMDEC);
946 dd_partition_system(fplog, step, cr,
947 bMasterState, nstglobalcomm,
948 state_global, top_global, ir,
949 state, &f, mdatoms, top, fr,
950 vsite, shellfc, constr,
952 do_verbose && !bPMETuneRunning);
953 wallcycle_stop(wcycle, ewcDOMDEC);
954 /* If using an iterative integrator, reallocate space to match the decomposition */
958 if (MASTER(cr) && do_log)
960 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
963 if (ir->efep != efepNO)
965 update_mdatoms(mdatoms, state->lambda[efptMASS]);
968 if ((bRerunMD && rerun_fr.bV) || bExchanged)
971 /* We need the kinetic energy at minus the half step for determining
972 * the full step kinetic energy and possibly for T-coupling.*/
973 /* This may not be quite working correctly yet . . . . */
974 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
975 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
976 constr, NULL, FALSE, state->box,
977 top_global, &bSumEkinhOld,
978 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
980 clear_mat(force_vir);
982 /* We write a checkpoint at this MD step when:
983 * either at an NS step when we signalled through gs,
984 * or at the last step (but not when we do not want confout),
985 * but never at the first step or with rerun.
987 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
988 (bLastStep && (Flags & MD_CONFOUT))) &&
989 step > ir->init_step && !bRerunMD);
992 gs.set[eglsCHKPT] = 0;
995 /* Determine the energy and pressure:
996 * at nstcalcenergy steps and at energy output steps (set below).
998 if (EI_VV(ir->eI) && (!bInitStep))
1000 /* for vv, the first half of the integration actually corresponds
1001 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1002 but the virial needs to be calculated on both the current step and the 'next' step. Future
1003 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1005 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1006 bCalcVir = bCalcEner ||
1007 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1011 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1012 bCalcVir = bCalcEner ||
1013 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1016 /* Do we need global communication ? */
1017 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1018 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1019 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1021 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1023 if (do_ene || do_log || bDoReplEx)
1030 /* these CGLO_ options remain the same throughout the iteration */
1031 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1032 (bGStat ? CGLO_GSTAT : 0)
1035 force_flags = (GMX_FORCE_STATECHANGED |
1036 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1037 GMX_FORCE_ALLFORCES |
1039 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1040 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1041 (bDoFEP ? GMX_FORCE_DHDL : 0)
1046 if (do_per_step(step, ir->nstcalclr))
1048 force_flags |= GMX_FORCE_DO_LR;
1054 /* Now is the time to relax the shells */
1055 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1056 ir, bNS, force_flags,
1059 state, f, force_vir, mdatoms,
1060 nrnb, wcycle, graph, groups,
1061 shellfc, fr, bBornRadii, t, mu_tot,
1063 mdoutf_get_fp_field(outf));
1073 /* The coordinates (x) are shifted (to get whole molecules)
1075 * This is parallellized as well, and does communication too.
1076 * Check comments in sim_util.c
1078 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1079 state->box, state->x, &state->hist,
1080 f, force_vir, mdatoms, enerd, fcd,
1081 state->lambda, graph,
1082 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1083 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1086 if (bVV && !bStartingFromCpt && !bRerunMD)
1087 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1091 wallcycle_start(wcycle, ewcUPDATE);
1092 if (ir->eI == eiVV && bInitStep)
1094 /* if using velocity verlet with full time step Ekin,
1095 * take the first half step only to compute the
1096 * virial for the first step. From there,
1097 * revert back to the initial coordinates
1098 * so that the input is actually the initial step.
1100 snew(vbuf, state->natoms);
1101 copy_rvecn(state->v, vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1105 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1106 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1109 /* If we are using twin-range interactions where the long-range component
1110 * is only evaluated every nstcalclr>1 steps, we should do a special update
1111 * step to combine the long-range forces on these steps.
1112 * For nstcalclr=1 this is not done, since the forces would have been added
1113 * directly to the short-range forces already.
1115 * TODO Remove various aspects of VV+twin-range in master
1116 * branch, because VV integrators did not ever support
1117 * twin-range multiple time stepping with constraints.
1119 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1121 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1122 f, bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1123 ekind, M, upd, bInitStep, etrtVELOCITY1,
1124 cr, nrnb, constr, &top->idef);
1126 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1128 gmx_iterate_init(&iterate, TRUE);
1130 /* for iterations, we save these vectors, as we will be self-consistently iterating
1133 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1135 /* save the state */
1136 if (iterate.bIterationActive)
1138 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1141 bFirstIterate = TRUE;
1142 while (bFirstIterate || iterate.bIterationActive)
1144 if (iterate.bIterationActive)
1146 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1147 if (bFirstIterate && bTrotter)
1149 /* The first time through, we need a decent first estimate
1150 of veta(t+dt) to compute the constraints. Do
1151 this by computing the box volume part of the
1152 trotter integration at this time. Nothing else
1153 should be changed by this routine here. If
1154 !(first time), we start with the previous value
1157 veta_save = state->veta;
1158 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1159 vetanew = state->veta;
1160 state->veta = veta_save;
1165 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1167 wallcycle_stop(wcycle, ewcUPDATE);
1168 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1169 state, fr->bMolPBC, graph, f,
1170 &top->idef, shake_vir,
1171 cr, nrnb, wcycle, upd, constr,
1172 TRUE, bCalcVir, vetanew);
1173 wallcycle_start(wcycle, ewcUPDATE);
1175 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1177 /* Correct the virial for multiple time stepping */
1178 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
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 wallcycle_stop(wcycle, ewcUPDATE);
1214 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1215 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1216 constr, NULL, FALSE, state->box,
1217 top_global, &bSumEkinhOld,
1220 | (bTemp ? CGLO_TEMPERATURE : 0)
1221 | (bPres ? CGLO_PRESSURE : 0)
1222 | (bPres ? CGLO_CONSTRAINT : 0)
1223 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1224 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1227 /* explanation of above:
1228 a) We compute Ekin at the full time step
1229 if 1) we are using the AveVel Ekin, and it's not the
1230 initial step, or 2) if we are using AveEkin, but need the full
1231 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1232 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1233 EkinAveVel because it's needed for the pressure */
1234 wallcycle_start(wcycle, ewcUPDATE);
1236 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1241 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1242 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1248 wallcycle_stop(wcycle, ewcUPDATE);
1249 /* We need the kinetic energy at minus the half step for determining
1250 * the full step kinetic energy and possibly for T-coupling.*/
1251 /* This may not be quite working correctly yet . . . . */
1252 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1253 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1254 constr, NULL, FALSE, state->box,
1255 top_global, &bSumEkinhOld,
1256 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1257 wallcycle_start(wcycle, ewcUPDATE);
1262 if (iterate.bIterationActive &&
1263 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1264 state->veta, &vetanew))
1268 bFirstIterate = FALSE;
1271 if (bTrotter && !bInitStep)
1273 copy_mat(shake_vir, state->svir_prev);
1274 copy_mat(force_vir, state->fvir_prev);
1275 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1277 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1278 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1279 enerd->term[F_EKIN] = trace(ekind->ekin);
1282 /* if it's the initial step, we performed this first step just to get the constraint virial */
1283 if (ir->eI == eiVV && bInitStep)
1285 copy_rvecn(vbuf, state->v, 0, state->natoms);
1288 wallcycle_stop(wcycle, ewcUPDATE);
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, 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_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1329 ir, state, state_global, top_global, fr,
1330 outf, mdebin, ekind, f, f_global,
1332 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1334 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1335 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1337 /* kludge -- virial is lost with restart for NPT control. Must restart */
1338 if (bStartingFromCpt && bVV)
1340 copy_mat(state->svir_prev, shake_vir);
1341 copy_mat(state->fvir_prev, force_vir);
1344 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1346 /* Check whether everything is still allright */
1347 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1348 #ifdef GMX_THREAD_MPI
1353 /* this is just make gs.sig compatible with the hack
1354 of sending signals around by MPI_Reduce with together with
1356 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1358 gs.sig[eglsSTOPCOND] = 1;
1360 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1362 gs.sig[eglsSTOPCOND] = -1;
1364 /* < 0 means stop at next step, > 0 means stop at next NS step */
1368 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1369 gmx_get_signal_name(),
1370 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1374 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1375 gmx_get_signal_name(),
1376 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1378 handled_stop_condition = (int)gmx_get_stop_condition();
1380 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1381 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1382 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1384 /* Signal to terminate the run */
1385 gs.sig[eglsSTOPCOND] = 1;
1388 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1390 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1393 if (bResetCountersHalfMaxH && MASTER(cr) &&
1394 elapsed_time > max_hours*60.0*60.0*0.495)
1396 gs.sig[eglsRESETCOUNTERS] = 1;
1399 if (ir->nstlist == -1 && !bRerunMD)
1401 /* When bGStatEveryStep=FALSE, global_stat is only called
1402 * when we check the atom displacements, not at NS steps.
1403 * This means that also the bonded interaction count check is not
1404 * performed immediately after NS. Therefore a few MD steps could
1405 * be performed with missing interactions.
1406 * But wrong energies are never written to file,
1407 * since energies are only written after global_stat
1410 if (step >= nlh.step_nscheck)
1412 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1413 nlh.scale_tot, state->x);
1417 /* This is not necessarily true,
1418 * but step_nscheck is determined quite conservatively.
1424 /* In parallel we only have to check for checkpointing in steps
1425 * where we do global communication,
1426 * otherwise the other nodes don't know.
1428 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1431 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1432 gs.set[eglsCHKPT] == 0)
1434 gs.sig[eglsCHKPT] = 1;
1437 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1442 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1444 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1446 gmx_bool bIfRandomize;
1447 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1448 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1449 if (constr && bIfRandomize)
1451 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1452 state, fr->bMolPBC, graph, f,
1453 &top->idef, tmp_vir,
1454 cr, nrnb, wcycle, upd, constr,
1455 TRUE, bCalcVir, vetanew);
1460 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1462 gmx_iterate_init(&iterate, TRUE);
1463 /* for iterations, we save these vectors, as we will be redoing the calculations */
1464 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1467 bFirstIterate = TRUE;
1468 while (bFirstIterate || iterate.bIterationActive)
1470 /* We now restore these vectors to redo the calculation with improved extended variables */
1471 if (iterate.bIterationActive)
1473 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1476 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1477 so scroll down for that logic */
1479 /* ######### START SECOND UPDATE STEP ################# */
1480 /* Box is changed in update() when we do pressure coupling,
1481 * but we should still use the old box for energy corrections and when
1482 * writing it to the energy file, so it matches the trajectory files for
1483 * the same timestep above. Make a copy in a separate array.
1485 copy_mat(state->box, lastbox);
1490 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1492 wallcycle_start(wcycle, ewcUPDATE);
1493 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1496 if (iterate.bIterationActive)
1504 /* we use a new value of scalevir to converge the iterations faster */
1505 scalevir = tracevir/trace(shake_vir);
1507 msmul(shake_vir, scalevir, shake_vir);
1508 m_add(force_vir, shake_vir, total_vir);
1509 clear_mat(shake_vir);
1511 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1512 /* We can only do Berendsen coupling after we have summed
1513 * the kinetic energy or virial. Since the happens
1514 * in global_state after update, we should only do it at
1515 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1520 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1521 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1526 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1528 /* velocity half-step update */
1529 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1530 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1531 ekind, M, upd, FALSE, etrtVELOCITY2,
1532 cr, nrnb, constr, &top->idef);
1535 /* Above, initialize just copies ekinh into ekin,
1536 * it doesn't copy position (for VV),
1537 * and entire integrator for MD.
1540 if (ir->eI == eiVVAK)
1542 /* We probably only need md->homenr, not state->natoms */
1543 if (state->natoms > cbuf_nalloc)
1545 cbuf_nalloc = state->natoms;
1546 srenew(cbuf, cbuf_nalloc);
1548 copy_rvecn(state->x, cbuf, 0, state->natoms);
1550 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1552 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1553 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1554 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1555 wallcycle_stop(wcycle, ewcUPDATE);
1557 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1558 fr->bMolPBC, graph, f,
1559 &top->idef, shake_vir,
1560 cr, nrnb, wcycle, upd, constr,
1561 FALSE, bCalcVir, state->veta);
1563 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1565 /* Correct the virial for multiple time stepping */
1566 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1569 if (ir->eI == eiVVAK)
1571 /* erase F_EKIN and F_TEMP here? */
1572 /* just compute the kinetic energy at the half step to perform a trotter step */
1573 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1574 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1575 constr, NULL, FALSE, lastbox,
1576 top_global, &bSumEkinhOld,
1577 cglo_flags | CGLO_TEMPERATURE
1579 wallcycle_start(wcycle, ewcUPDATE);
1580 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1581 /* now we know the scaling, we can compute the positions again again */
1582 copy_rvecn(cbuf, state->x, 0, state->natoms);
1584 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1586 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1587 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1588 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1589 wallcycle_stop(wcycle, ewcUPDATE);
1591 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1592 /* are the small terms in the shake_vir here due
1593 * to numerical errors, or are they important
1594 * physically? I'm thinking they are just errors, but not completely sure.
1595 * For now, will call without actually constraining, constr=NULL*/
1596 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1597 state, fr->bMolPBC, graph, f,
1598 &top->idef, tmp_vir,
1599 cr, nrnb, wcycle, upd, NULL,
1605 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1608 if (fr->bSepDVDL && fplog && do_log)
1610 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1614 /* this factor or 2 correction is necessary
1615 because half of the constraint force is removed
1616 in the vv step, so we have to double it. See
1617 the Redmine issue #1255. It is not yet clear
1618 if the factor of 2 is exact, or just a very
1619 good approximation, and this will be
1620 investigated. The next step is to see if this
1621 can be done adding a dhdl contribution from the
1622 rattle step, but this is somewhat more
1623 complicated with the current code. Will be
1624 investigated, hopefully for 4.6.3. However,
1625 this current solution is much better than
1626 having it completely wrong.
1628 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1632 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1637 /* Need to unshift here */
1638 unshift_self(graph, state->box, state->x);
1643 wallcycle_start(wcycle, ewcVSITECONSTR);
1646 shift_self(graph, state->box, state->x);
1648 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1649 top->idef.iparams, top->idef.il,
1650 fr->ePBC, fr->bMolPBC, cr, state->box);
1654 unshift_self(graph, state->box, state->x);
1656 wallcycle_stop(wcycle, ewcVSITECONSTR);
1659 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1660 /* With Leap-Frog we can skip compute_globals at
1661 * non-communication steps, but we need to calculate
1662 * the kinetic energy one step before communication.
1664 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1666 if (ir->nstlist == -1 && bFirstIterate)
1668 gs.sig[eglsNABNSB] = nlh.nabnsb;
1670 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1671 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1673 bFirstIterate ? &gs : NULL,
1674 (step_rel % gs.nstms == 0) &&
1675 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1677 top_global, &bSumEkinhOld,
1679 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1680 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1681 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1682 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1683 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1684 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1687 if (ir->nstlist == -1 && bFirstIterate)
1689 nlh.nabnsb = gs.set[eglsNABNSB];
1690 gs.set[eglsNABNSB] = 0;
1693 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1694 /* ############# END CALC EKIN AND PRESSURE ################# */
1696 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1697 the virial that should probably be addressed eventually. state->veta has better properies,
1698 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1699 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1701 if (iterate.bIterationActive &&
1702 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1703 trace(shake_vir), &tracevir))
1707 bFirstIterate = FALSE;
1710 if (!bVV || bRerunMD)
1712 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1713 sum_dhdl(enerd, state->lambda, ir->fepvals);
1715 update_box(fplog, step, ir, mdatoms, state, f,
1716 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1718 /* ################# END UPDATE STEP 2 ################# */
1719 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1721 /* The coordinates (x) were unshifted in update */
1724 /* We will not sum ekinh_old,
1725 * so signal that we still have to do it.
1727 bSumEkinhOld = TRUE;
1730 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1732 /* use the directly determined last velocity, not actually the averaged half steps */
1733 if (bTrotter && ir->eI == eiVV)
1735 enerd->term[F_EKIN] = last_ekin;
1737 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1741 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1745 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1747 /* ######### END PREPARING EDR OUTPUT ########### */
1752 gmx_bool do_dr, do_or;
1754 if (fplog && do_log && bDoExpanded)
1756 /* only needed if doing expanded ensemble */
1757 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1758 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1760 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1764 upd_mdebin(mdebin, bDoDHDL, TRUE,
1765 t, mdatoms->tmass, enerd, state,
1766 ir->fepvals, ir->expandedvals, lastbox,
1767 shake_vir, force_vir, total_vir, pres,
1768 ekind, mu_tot, constr);
1772 upd_mdebin_step(mdebin);
1775 do_dr = do_per_step(step, ir->nstdisreout);
1776 do_or = do_per_step(step, ir->nstorireout);
1778 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1780 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1782 if (ir->ePull != epullNO)
1784 pull_print_output(ir->pull, step, t);
1787 if (do_per_step(step, ir->nstlog))
1789 if (fflush(fplog) != 0)
1791 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1797 /* Have to do this part _after_ outputting the logfile and the edr file */
1798 /* Gets written into the state at the beginning of next loop*/
1799 state->fep_state = lamnew;
1801 /* Print the remaining wall clock time for the run */
1802 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1806 fprintf(stderr, "\n");
1808 print_time(stderr, walltime_accounting, step, ir, cr);
1811 /* Ion/water position swapping.
1812 * Not done in last step since trajectory writing happens before this call
1813 * in the MD loop and exchanges would be lost anyway. */
1814 bNeedRepartition = FALSE;
1815 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1816 do_per_step(step, ir->swap->nstswap))
1818 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1819 bRerunMD ? rerun_fr.x : state->x,
1820 bRerunMD ? rerun_fr.box : state->box,
1821 top_global, MASTER(cr) && bVerbose, bRerunMD);
1823 if (bNeedRepartition && DOMAINDECOMP(cr))
1825 dd_collect_state(cr->dd, state, state_global);
1829 /* Replica exchange */
1833 bExchanged = replica_exchange(fplog, cr, repl_ex,
1834 state_global, enerd,
1838 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1840 dd_partition_system(fplog, step, cr, TRUE, 1,
1841 state_global, top_global, ir,
1842 state, &f, mdatoms, top, fr,
1843 vsite, shellfc, constr,
1844 nrnb, wcycle, FALSE);
1849 bStartingFromCpt = FALSE;
1851 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1852 /* With all integrators, except VV, we need to retain the pressure
1853 * at the current step for coupling at the next step.
1855 if ((state->flags & (1<<estPRES_PREV)) &&
1857 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1859 /* Store the pressure in t_state for pressure coupling
1860 * at the next MD step.
1862 copy_mat(pres, state->pres_prev);
1865 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1867 if ( (membed != NULL) && (!bLastStep) )
1869 rescale_membed(step_rel, membed, state_global->x);
1876 /* read next frame from input trajectory */
1877 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1882 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1886 if (!bRerunMD || !rerun_fr.bStep)
1888 /* increase the MD step number */
1893 cycles = wallcycle_stop(wcycle, ewcSTEP);
1894 if (DOMAINDECOMP(cr) && wcycle)
1896 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1899 if (bPMETuneRunning || bPMETuneTry)
1901 /* PME grid + cut-off optimization with GPUs or PME nodes */
1903 /* Count the total cycles over the last steps */
1904 cycles_pmes += cycles;
1906 /* We can only switch cut-off at NS steps */
1907 if (step % ir->nstlist == 0)
1909 /* PME grid + cut-off optimization with GPUs or PME nodes */
1912 if (DDMASTER(cr->dd))
1914 /* PME node load is too high, start tuning */
1915 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1917 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1919 if (bPMETuneRunning &&
1920 fr->nbv->bUseGPU && DOMAINDECOMP(cr) &&
1921 !(cr->duty & DUTY_PME))
1923 /* Lock DLB=auto to off (does nothing when DLB=yes/no).
1924 * With GPUs + separate PME ranks, we don't want DLB.
1925 * This could happen when we scan coarse grids and
1926 * it would then never be turned off again.
1927 * This would hurt performance at the final, optimal
1928 * grid spacing, where DLB almost never helps.
1929 * Also, DLB can limit the cut-off for PME tuning.
1931 dd_dlb_set_lock(cr->dd, TRUE);
1934 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1936 bPMETuneTry = FALSE;
1939 if (bPMETuneRunning)
1941 /* init_step might not be a multiple of nstlist,
1942 * but the first cycle is always skipped anyhow.
1945 pme_load_balance(pme_loadbal, cr,
1946 (bVerbose && MASTER(cr)) ? stderr : NULL,
1948 ir, state, cycles_pmes,
1949 fr->ic, fr->nbv, &fr->pmedata,
1952 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1953 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1954 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1955 fr->rlist = fr->ic->rlist;
1956 fr->rlistlong = fr->ic->rlistlong;
1957 fr->rcoulomb = fr->ic->rcoulomb;
1958 fr->rvdw = fr->ic->rvdw;
1960 if (ir->eDispCorr != edispcNO)
1962 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1965 if (!bPMETuneRunning &&
1967 dd_dlb_is_locked(cr->dd))
1969 /* Unlock the DLB=auto, DLB is allowed to activate
1970 * (but we don't expect it to activate in most cases).
1972 dd_dlb_set_lock(cr->dd, FALSE);
1979 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1980 gs.set[eglsRESETCOUNTERS] != 0)
1982 /* Reset all the counters related to performance over the run */
1983 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1984 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1985 wcycle_set_reset_counters(wcycle, -1);
1986 if (!(cr->duty & DUTY_PME))
1988 /* Tell our PME node to reset its counters */
1989 gmx_pme_send_resetcounters(cr, step);
1991 /* Correct max_hours for the elapsed time */
1992 max_hours -= elapsed_time/(60.0*60.0);
1993 bResetCountersHalfMaxH = FALSE;
1994 gs.set[eglsRESETCOUNTERS] = 0;
1997 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1998 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
2001 /* End of main MD loop */
2004 /* Closing TNG files can include compressing data. Therefore it is good to do that
2005 * before stopping the time measurements. */
2006 mdoutf_tng_close(outf);
2008 /* Stop measuring walltime */
2009 walltime_accounting_end(walltime_accounting);
2011 if (bRerunMD && MASTER(cr))
2016 if (!(cr->duty & DUTY_PME))
2018 /* Tell the PME only node to finish */
2019 gmx_pme_send_finish(cr);
2024 if (ir->nstcalcenergy > 0 && !bRerunMD)
2026 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
2027 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
2034 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2036 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)));
2037 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
2040 if (pme_loadbal != NULL)
2042 pme_loadbal_done(pme_loadbal, cr, fplog,
2043 fr->nbv != NULL && fr->nbv->bUseGPU);
2046 if (shellfc && fplog)
2048 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
2049 (nconverged*100.0)/step_rel);
2050 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
2054 if (repl_ex_nst > 0 && MASTER(cr))
2056 print_replica_exchange_statistics(fplog, repl_ex);
2059 /* IMD cleanup, if bIMD is TRUE. */
2060 IMD_finalize(ir->bIMD, ir->imd);
2062 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);