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.
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;
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 gmx_bool bIMDstep = FALSE;
228 /* Temporary addition for FAHCORE checkpointing */
232 /* Check for special mdrun options */
233 bRerunMD = (Flags & MD_RERUN);
234 bAppend = (Flags & MD_APPENDFILES);
235 if (Flags & MD_RESETCOUNTERSHALFWAY)
239 /* Signal to reset the counters half the simulation steps. */
240 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
242 /* Signal to reset the counters halfway the simulation time. */
243 bResetCountersHalfMaxH = (max_hours > 0);
246 /* md-vv uses averaged full step velocities for T-control
247 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
248 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
250 if (bVV) /* to store the initial velocities while computing virial */
252 snew(cbuf, top_global->natoms);
254 /* all the iteratative cases - only if there are constraints */
255 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
256 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
257 false in this step. The correct value, true or false,
258 is set at each step, as it depends on the frequency of temperature
259 and pressure control.*/
260 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
264 /* Since we don't know if the frames read are related in any way,
265 * rebuild the neighborlist at every step.
268 ir->nstcalcenergy = 1;
272 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
274 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
275 bGStatEveryStep = (nstglobalcomm == 1);
277 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
280 "To reduce the energy communication with nstlist = -1\n"
281 "the neighbor list validity should not be checked at every step,\n"
282 "this means that exact integration is not guaranteed.\n"
283 "The neighbor list validity is checked after:\n"
284 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
285 "In most cases this will result in exact integration.\n"
286 "This reduces the energy communication by a factor of 2 to 3.\n"
287 "If you want less energy communication, set nstlist > 3.\n\n");
292 ir->nstxout_compressed = 0;
294 groups = &top_global->groups;
297 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
298 &(state_global->fep_state), lam0,
299 nrnb, top_global, &upd,
300 nfile, fnm, &outf, &mdebin,
301 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags, wcycle);
303 clear_mat(total_vir);
305 /* Energy terms and groups */
307 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
309 if (DOMAINDECOMP(cr))
315 snew(f, top_global->natoms);
318 /* Kinetic energy data */
320 init_ekindata(fplog, top_global, &(ir->opts), ekind);
321 /* needed for iteration of constraints */
323 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
324 /* Copy the cos acceleration to the groups struct */
325 ekind->cosacc.cos_accel = ir->cos_accel;
327 gstat = global_stat_init(ir);
330 /* Check for polarizable models and flexible constraints */
331 shellfc = init_shell_flexcon(fplog,
332 top_global, n_flexible_constraints(constr),
333 (ir->bContinuation ||
334 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
335 NULL : state_global->x);
336 if (shellfc && ir->nstcalcenergy != 1)
338 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);
340 if (shellfc && DOMAINDECOMP(cr))
342 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
344 if (shellfc && ir->eI == eiNM)
346 /* Currently shells don't work with Normal Modes */
347 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");
350 if (vsite && ir->eI == eiNM)
352 /* Currently virtual sites don't work with Normal Modes */
353 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");
358 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
359 set_deform_reference_box(upd,
360 deform_init_init_step_tpx,
361 deform_init_box_tpx);
362 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
366 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
367 if ((io > 2000) && MASTER(cr))
370 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
375 if (DOMAINDECOMP(cr))
377 top = dd_init_local_top(top_global);
380 dd_init_local_state(cr->dd, state_global, state);
382 if (DDMASTER(cr->dd) && ir->nstfout)
384 snew(f_global, state_global->natoms);
389 top = gmx_mtop_generate_local_top(top_global, ir);
391 forcerec_set_excl_load(fr, top);
393 state = serial_init_local_state(state_global);
396 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
400 set_vsite_top(vsite, top, mdatoms, cr);
403 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
405 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
410 make_local_shells(cr, mdatoms, shellfc);
413 setup_bonded_threading(fr, &top->idef);
416 /* Set up interactive MD (IMD) */
417 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
418 nfile, fnm, oenv, imdport, Flags);
420 if (DOMAINDECOMP(cr))
422 /* Distribute the charge groups over the nodes from the master node */
423 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
424 state_global, top_global, ir,
425 state, &f, mdatoms, top, fr,
426 vsite, shellfc, constr,
427 nrnb, wcycle, FALSE);
431 update_mdatoms(mdatoms, state->lambda[efptMASS]);
433 if (opt2bSet("-cpi", nfile, fnm))
435 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
439 bStateFromCP = FALSE;
444 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
451 /* Update mdebin with energy history if appending to output files */
452 if (Flags & MD_APPENDFILES)
454 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
458 /* We might have read an energy history from checkpoint,
459 * free the allocated memory and reset the counts.
461 done_energyhistory(&state_global->enerhist);
462 init_energyhistory(&state_global->enerhist);
465 /* Set the initial energy history in state by updating once */
466 update_energyhistory(&state_global->enerhist, mdebin);
469 /* Initialize constraints */
470 if (constr && !DOMAINDECOMP(cr))
472 set_constraints(constr, top, ir, mdatoms, cr);
475 if (repl_ex_nst > 0 && MASTER(cr))
477 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
478 repl_ex_nst, repl_ex_nex, repl_ex_seed);
481 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
482 * PME tuning is not supported with PME only for LJ and not for Coulomb.
484 if ((Flags & MD_TUNEPME) &&
485 EEL_PME(fr->eeltype) &&
486 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
489 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
491 if (cr->duty & DUTY_PME)
493 /* Start tuning right away, as we can't measure the load */
494 bPMETuneRunning = TRUE;
498 /* Separate PME nodes, we can measure the PP/PME load balance */
503 if (!ir->bContinuation && !bRerunMD)
505 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
507 /* Set the velocities of frozen particles to zero */
508 for (i = 0; i < mdatoms->homenr; i++)
510 for (m = 0; m < DIM; m++)
512 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
522 /* Constrain the initial coordinates and velocities */
523 do_constrain_first(fplog, constr, ir, mdatoms, state,
528 /* Construct the virtual sites for the initial configuration */
529 construct_vsites(vsite, state->x, ir->delta_t, NULL,
530 top->idef.iparams, top->idef.il,
531 fr->ePBC, fr->bMolPBC, cr, state->box);
537 /* set free energy calculation frequency as the minimum
538 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
539 nstfep = ir->fepvals->nstdhdl;
542 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
546 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
549 /* I'm assuming we need global communication the first time! MRS */
550 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
551 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
552 | (bVV ? CGLO_PRESSURE : 0)
553 | (bVV ? CGLO_CONSTRAINT : 0)
554 | (bRerunMD ? CGLO_RERUNMD : 0)
555 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
557 bSumEkinhOld = FALSE;
558 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
559 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
560 constr, NULL, FALSE, state->box,
561 top_global, &bSumEkinhOld, cglo_flags);
562 if (ir->eI == eiVVAK)
564 /* a second call to get the half step temperature initialized as well */
565 /* we do the same call as above, but turn the pressure off -- internally to
566 compute_globals, this is recognized as a velocity verlet half-step
567 kinetic energy calculation. This minimized excess variables, but
568 perhaps loses some logic?*/
570 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
571 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
572 constr, NULL, FALSE, state->box,
573 top_global, &bSumEkinhOld,
574 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
577 /* Calculate the initial half step temperature, and save the ekinh_old */
578 if (!(Flags & MD_STARTFROMCPT))
580 for (i = 0; (i < ir->opts.ngtc); i++)
582 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
587 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
588 and there is no previous step */
591 /* if using an iterative algorithm, we need to create a working directory for the state. */
594 bufstate = init_bufstate(state);
597 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
598 temperature control */
599 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
603 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
606 "RMS relative constraint deviation after constraining: %.2e\n",
607 constr_rmsd(constr, FALSE));
609 if (EI_STATE_VELOCITY(ir->eI))
611 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
615 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
616 " input trajectory '%s'\n\n",
617 *(top_global->name), opt2fn("-rerun", nfile, fnm));
620 fprintf(stderr, "Calculated time to finish depends on nsteps from "
621 "run input file,\nwhich may not correspond to the time "
622 "needed to process input trajectory.\n\n");
628 fprintf(stderr, "starting mdrun '%s'\n",
629 *(top_global->name));
632 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
636 sprintf(tbuf, "%s", "infinite");
638 if (ir->init_step > 0)
640 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
641 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
642 gmx_step_str(ir->init_step, sbuf2),
643 ir->init_step*ir->delta_t);
647 fprintf(stderr, "%s steps, %s ps.\n",
648 gmx_step_str(ir->nsteps, sbuf), tbuf);
651 fprintf(fplog, "\n");
654 walltime_accounting_start(walltime_accounting);
655 wallcycle_start(wcycle, ewcRUN);
656 print_start(fplog, cr, walltime_accounting, "mdrun");
658 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
660 chkpt_ret = fcCheckPointParallel( cr->nodeid,
664 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
669 /***********************************************************
673 ************************************************************/
675 /* if rerunMD then read coordinates and velocities from input trajectory */
678 if (getenv("GMX_FORCE_UPDATE"))
686 bNotLastFrame = read_first_frame(oenv, &status,
687 opt2fn("-rerun", nfile, fnm),
688 &rerun_fr, TRX_NEED_X | TRX_READ_V);
689 if (rerun_fr.natoms != top_global->natoms)
692 "Number of atoms in trajectory (%d) does not match the "
693 "run input file (%d)\n",
694 rerun_fr.natoms, top_global->natoms);
696 if (ir->ePBC != epbcNONE)
700 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);
702 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
704 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
711 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
714 if (ir->ePBC != epbcNONE)
716 /* Set the shift vectors.
717 * Necessary here when have a static box different from the tpr box.
719 calc_shifts(rerun_fr.box, fr->shift_vec);
723 /* loop over MD steps or if rerunMD to end of input trajectory */
725 /* Skip the first Nose-Hoover integration when we get the state from tpx */
726 bStateFromTPX = !bStateFromCP;
727 bInitStep = bFirstStep && (bStateFromTPX || bVV);
728 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
730 bSumEkinhOld = FALSE;
733 bNeedRepartition = FALSE;
735 init_global_signals(&gs, cr, ir, repl_ex_nst);
737 step = ir->init_step;
740 if (ir->nstlist == -1)
742 init_nlistheuristics(&nlh, bGStatEveryStep, step);
745 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
747 /* check how many steps are left in other sims */
748 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
752 /* and stop now if we should */
753 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
754 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
755 while (!bLastStep || (bRerunMD && bNotLastFrame))
758 wallcycle_start(wcycle, ewcSTEP);
764 step = rerun_fr.step;
765 step_rel = step - ir->init_step;
778 bLastStep = (step_rel == ir->nsteps);
779 t = t0 + step*ir->delta_t;
782 if (ir->efep != efepNO || ir->bSimTemp)
784 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
785 requiring different logic. */
787 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
788 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
789 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
790 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
791 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
794 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
795 do_per_step(step, repl_ex_nst));
799 update_annealing_target_temp(&(ir->opts), t);
804 if (!DOMAINDECOMP(cr) || MASTER(cr))
806 for (i = 0; i < state_global->natoms; i++)
808 copy_rvec(rerun_fr.x[i], state_global->x[i]);
812 for (i = 0; i < state_global->natoms; i++)
814 copy_rvec(rerun_fr.v[i], state_global->v[i]);
819 for (i = 0; i < state_global->natoms; i++)
821 clear_rvec(state_global->v[i]);
825 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
826 " Ekin, temperature and pressure are incorrect,\n"
827 " the virial will be incorrect when constraints are present.\n"
829 bRerunWarnNoV = FALSE;
833 copy_mat(rerun_fr.box, state_global->box);
834 copy_mat(state_global->box, state->box);
836 if (vsite && (Flags & MD_RERUN_VSITE))
838 if (DOMAINDECOMP(cr))
840 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
844 /* Following is necessary because the graph may get out of sync
845 * with the coordinates if we only have every N'th coordinate set
847 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
848 shift_self(graph, state->box, state->x);
850 construct_vsites(vsite, state->x, ir->delta_t, state->v,
851 top->idef.iparams, top->idef.il,
852 fr->ePBC, fr->bMolPBC, cr, state->box);
855 unshift_self(graph, state->box, state->x);
860 /* Stop Center of Mass motion */
861 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
865 /* for rerun MD always do Neighbour Searching */
866 bNS = (bFirstStep || ir->nstlist != 0);
871 /* Determine whether or not to do Neighbour Searching and LR */
872 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
874 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
875 (ir->nstlist == -1 && nlh.nabnsb > 0));
877 if (bNS && ir->nstlist == -1)
879 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
883 /* check whether we should stop because another simulation has
887 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
888 (multisim_nsteps != ir->nsteps) )
895 "Stopping simulation %d because another one has finished\n",
899 gs.sig[eglsCHKPT] = 1;
904 /* < 0 means stop at next step, > 0 means stop at next NS step */
905 if ( (gs.set[eglsSTOPCOND] < 0) ||
906 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
911 /* Determine whether or not to update the Born radii if doing GB */
912 bBornRadii = bFirstStep;
913 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
918 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
919 do_verbose = bVerbose &&
920 (step % stepout == 0 || bFirstStep || bLastStep);
922 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
930 bMasterState = FALSE;
931 /* Correct the new box if it is too skewed */
932 if (DYNAMIC_BOX(*ir))
934 if (correct_box(fplog, step, state->box, graph))
939 if (DOMAINDECOMP(cr) && bMasterState)
941 dd_collect_state(cr->dd, state, state_global);
945 if (DOMAINDECOMP(cr))
947 /* Repartition the domain decomposition */
948 wallcycle_start(wcycle, ewcDOMDEC);
949 dd_partition_system(fplog, step, cr,
950 bMasterState, nstglobalcomm,
951 state_global, top_global, ir,
952 state, &f, mdatoms, top, fr,
953 vsite, shellfc, constr,
955 do_verbose && !bPMETuneRunning);
956 wallcycle_stop(wcycle, ewcDOMDEC);
957 /* If using an iterative integrator, reallocate space to match the decomposition */
961 if (MASTER(cr) && do_log)
963 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
966 if (ir->efep != efepNO)
968 update_mdatoms(mdatoms, state->lambda[efptMASS]);
971 if ((bRerunMD && rerun_fr.bV) || bExchanged)
974 /* We need the kinetic energy at minus the half step for determining
975 * the full step kinetic energy and possibly for T-coupling.*/
976 /* This may not be quite working correctly yet . . . . */
977 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
978 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
979 constr, NULL, FALSE, state->box,
980 top_global, &bSumEkinhOld,
981 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
983 clear_mat(force_vir);
985 /* We write a checkpoint at this MD step when:
986 * either at an NS step when we signalled through gs,
987 * or at the last step (but not when we do not want confout),
988 * but never at the first step or with rerun.
990 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
991 (bLastStep && (Flags & MD_CONFOUT))) &&
992 step > ir->init_step && !bRerunMD);
995 gs.set[eglsCHKPT] = 0;
998 /* Determine the energy and pressure:
999 * at nstcalcenergy steps and at energy output steps (set below).
1001 if (EI_VV(ir->eI) && (!bInitStep))
1003 /* for vv, the first half of the integration actually corresponds
1004 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1005 but the virial needs to be calculated on both the current step and the 'next' step. Future
1006 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1008 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1009 bCalcVir = bCalcEner ||
1010 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1014 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1015 bCalcVir = bCalcEner ||
1016 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1019 /* Do we need global communication ? */
1020 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1021 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1022 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1024 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1026 if (do_ene || do_log || bDoReplEx)
1033 /* these CGLO_ options remain the same throughout the iteration */
1034 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1035 (bGStat ? CGLO_GSTAT : 0)
1038 force_flags = (GMX_FORCE_STATECHANGED |
1039 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1040 GMX_FORCE_ALLFORCES |
1042 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1043 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1044 (bDoFEP ? GMX_FORCE_DHDL : 0)
1049 if (do_per_step(step, ir->nstcalclr))
1051 force_flags |= GMX_FORCE_DO_LR;
1057 /* Now is the time to relax the shells */
1058 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1059 ir, bNS, force_flags,
1062 state, f, force_vir, mdatoms,
1063 nrnb, wcycle, graph, groups,
1064 shellfc, fr, bBornRadii, t, mu_tot,
1066 mdoutf_get_fp_field(outf));
1076 /* The coordinates (x) are shifted (to get whole molecules)
1078 * This is parallellized as well, and does communication too.
1079 * Check comments in sim_util.c
1081 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1082 state->box, state->x, &state->hist,
1083 f, force_vir, mdatoms, enerd, fcd,
1084 state->lambda, graph,
1085 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1086 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1089 if (bVV && !bStartingFromCpt && !bRerunMD)
1090 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1092 wallcycle_start(wcycle, ewcUPDATE);
1093 if (ir->eI == eiVV && bInitStep)
1095 /* if using velocity verlet with full time step Ekin,
1096 * take the first half step only to compute the
1097 * virial for the first step. From there,
1098 * revert back to the initial coordinates
1099 * so that the input is actually the initial step.
1101 copy_rvecn(state->v, cbuf, 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 (bInitStep && ir->eI == eiVV)
1285 copy_rvecn(cbuf, state->v, 0, state->natoms);
1287 wallcycle_stop(wcycle, ewcUPDATE);
1290 /* MRS -- now done iterating -- compute the conserved quantity */
1293 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1296 last_ekin = enerd->term[F_EKIN];
1298 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1300 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1302 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1305 sum_dhdl(enerd, state->lambda, ir->fepvals);
1309 /* ######## END FIRST UPDATE STEP ############## */
1310 /* ######## If doing VV, we now have v(dt) ###### */
1313 /* perform extended ensemble sampling in lambda - we don't
1314 actually move to the new state before outputting
1315 statistics, but if performing simulated tempering, we
1316 do update the velocities and the tau_t. */
1318 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1319 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1320 copy_df_history(&state_global->dfhist, &state->dfhist);
1323 /* Now we have the energies and forces corresponding to the
1324 * coordinates at time t. We must output all of this before
1327 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1328 ir, state, state_global, top_global, fr,
1329 outf, mdebin, ekind, f, f_global,
1331 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1333 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1334 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1336 /* kludge -- virial is lost with restart for NPT control. Must restart */
1337 if (bStartingFromCpt && bVV)
1339 copy_mat(state->svir_prev, shake_vir);
1340 copy_mat(state->fvir_prev, force_vir);
1343 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1345 /* Check whether everything is still allright */
1346 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1347 #ifdef GMX_THREAD_MPI
1352 /* this is just make gs.sig compatible with the hack
1353 of sending signals around by MPI_Reduce with together with
1355 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1357 gs.sig[eglsSTOPCOND] = 1;
1359 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1361 gs.sig[eglsSTOPCOND] = -1;
1363 /* < 0 means stop at next step, > 0 means stop at next NS step */
1367 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1368 gmx_get_signal_name(),
1369 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1373 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1374 gmx_get_signal_name(),
1375 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1377 handled_stop_condition = (int)gmx_get_stop_condition();
1379 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1380 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1381 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1383 /* Signal to terminate the run */
1384 gs.sig[eglsSTOPCOND] = 1;
1387 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1389 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1392 if (bResetCountersHalfMaxH && MASTER(cr) &&
1393 elapsed_time > max_hours*60.0*60.0*0.495)
1395 gs.sig[eglsRESETCOUNTERS] = 1;
1398 if (ir->nstlist == -1 && !bRerunMD)
1400 /* When bGStatEveryStep=FALSE, global_stat is only called
1401 * when we check the atom displacements, not at NS steps.
1402 * This means that also the bonded interaction count check is not
1403 * performed immediately after NS. Therefore a few MD steps could
1404 * be performed with missing interactions.
1405 * But wrong energies are never written to file,
1406 * since energies are only written after global_stat
1409 if (step >= nlh.step_nscheck)
1411 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1412 nlh.scale_tot, state->x);
1416 /* This is not necessarily true,
1417 * but step_nscheck is determined quite conservatively.
1423 /* In parallel we only have to check for checkpointing in steps
1424 * where we do global communication,
1425 * otherwise the other nodes don't know.
1427 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1430 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1431 gs.set[eglsCHKPT] == 0)
1433 gs.sig[eglsCHKPT] = 1;
1436 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1441 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1443 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1445 gmx_bool bIfRandomize;
1446 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1447 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1448 if (constr && bIfRandomize)
1450 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1451 state, fr->bMolPBC, graph, f,
1452 &top->idef, tmp_vir,
1453 cr, nrnb, wcycle, upd, constr,
1454 TRUE, bCalcVir, vetanew);
1459 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1461 gmx_iterate_init(&iterate, TRUE);
1462 /* for iterations, we save these vectors, as we will be redoing the calculations */
1463 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1466 bFirstIterate = TRUE;
1467 while (bFirstIterate || iterate.bIterationActive)
1469 /* We now restore these vectors to redo the calculation with improved extended variables */
1470 if (iterate.bIterationActive)
1472 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1475 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1476 so scroll down for that logic */
1478 /* ######### START SECOND UPDATE STEP ################# */
1479 /* Box is changed in update() when we do pressure coupling,
1480 * but we should still use the old box for energy corrections and when
1481 * writing it to the energy file, so it matches the trajectory files for
1482 * the same timestep above. Make a copy in a separate array.
1484 copy_mat(state->box, lastbox);
1489 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1491 wallcycle_start(wcycle, ewcUPDATE);
1492 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1495 if (iterate.bIterationActive)
1503 /* we use a new value of scalevir to converge the iterations faster */
1504 scalevir = tracevir/trace(shake_vir);
1506 msmul(shake_vir, scalevir, shake_vir);
1507 m_add(force_vir, shake_vir, total_vir);
1508 clear_mat(shake_vir);
1510 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1511 /* We can only do Berendsen coupling after we have summed
1512 * the kinetic energy or virial. Since the happens
1513 * in global_state after update, we should only do it at
1514 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1519 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1520 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1525 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1527 /* velocity half-step update */
1528 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1529 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1530 ekind, M, upd, FALSE, etrtVELOCITY2,
1531 cr, nrnb, constr, &top->idef);
1534 /* Above, initialize just copies ekinh into ekin,
1535 * it doesn't copy position (for VV),
1536 * and entire integrator for MD.
1539 if (ir->eI == eiVVAK)
1541 copy_rvecn(state->x, cbuf, 0, state->natoms);
1543 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1545 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1546 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1547 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1548 wallcycle_stop(wcycle, ewcUPDATE);
1550 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1551 fr->bMolPBC, graph, f,
1552 &top->idef, shake_vir,
1553 cr, nrnb, wcycle, upd, constr,
1554 FALSE, bCalcVir, state->veta);
1556 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1558 /* Correct the virial for multiple time stepping */
1559 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1562 if (ir->eI == eiVVAK)
1564 /* erase F_EKIN and F_TEMP here? */
1565 /* just compute the kinetic energy at the half step to perform a trotter step */
1566 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1567 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1568 constr, NULL, FALSE, lastbox,
1569 top_global, &bSumEkinhOld,
1570 cglo_flags | CGLO_TEMPERATURE
1572 wallcycle_start(wcycle, ewcUPDATE);
1573 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1574 /* now we know the scaling, we can compute the positions again again */
1575 copy_rvecn(cbuf, state->x, 0, state->natoms);
1577 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1579 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1580 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1581 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1582 wallcycle_stop(wcycle, ewcUPDATE);
1584 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1585 /* are the small terms in the shake_vir here due
1586 * to numerical errors, or are they important
1587 * physically? I'm thinking they are just errors, but not completely sure.
1588 * For now, will call without actually constraining, constr=NULL*/
1589 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1590 state, fr->bMolPBC, graph, f,
1591 &top->idef, tmp_vir,
1592 cr, nrnb, wcycle, upd, NULL,
1598 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1601 if (fr->bSepDVDL && fplog && do_log)
1603 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1607 /* this factor or 2 correction is necessary
1608 because half of the constraint force is removed
1609 in the vv step, so we have to double it. See
1610 the Redmine issue #1255. It is not yet clear
1611 if the factor of 2 is exact, or just a very
1612 good approximation, and this will be
1613 investigated. The next step is to see if this
1614 can be done adding a dhdl contribution from the
1615 rattle step, but this is somewhat more
1616 complicated with the current code. Will be
1617 investigated, hopefully for 4.6.3. However,
1618 this current solution is much better than
1619 having it completely wrong.
1621 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1625 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1630 /* Need to unshift here */
1631 unshift_self(graph, state->box, state->x);
1636 wallcycle_start(wcycle, ewcVSITECONSTR);
1639 shift_self(graph, state->box, state->x);
1641 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1642 top->idef.iparams, top->idef.il,
1643 fr->ePBC, fr->bMolPBC, cr, state->box);
1647 unshift_self(graph, state->box, state->x);
1649 wallcycle_stop(wcycle, ewcVSITECONSTR);
1652 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1653 /* With Leap-Frog we can skip compute_globals at
1654 * non-communication steps, but we need to calculate
1655 * the kinetic energy one step before communication.
1657 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1659 if (ir->nstlist == -1 && bFirstIterate)
1661 gs.sig[eglsNABNSB] = nlh.nabnsb;
1663 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1664 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1666 bFirstIterate ? &gs : NULL,
1667 (step_rel % gs.nstms == 0) &&
1668 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1670 top_global, &bSumEkinhOld,
1672 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1673 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1674 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1675 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1676 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1677 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1680 if (ir->nstlist == -1 && bFirstIterate)
1682 nlh.nabnsb = gs.set[eglsNABNSB];
1683 gs.set[eglsNABNSB] = 0;
1686 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1687 /* ############# END CALC EKIN AND PRESSURE ################# */
1689 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1690 the virial that should probably be addressed eventually. state->veta has better properies,
1691 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1692 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1694 if (iterate.bIterationActive &&
1695 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1696 trace(shake_vir), &tracevir))
1700 bFirstIterate = FALSE;
1703 if (!bVV || bRerunMD)
1705 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1706 sum_dhdl(enerd, state->lambda, ir->fepvals);
1708 update_box(fplog, step, ir, mdatoms, state, f,
1709 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1711 /* ################# END UPDATE STEP 2 ################# */
1712 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1714 /* The coordinates (x) were unshifted in update */
1717 /* We will not sum ekinh_old,
1718 * so signal that we still have to do it.
1720 bSumEkinhOld = TRUE;
1723 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1725 /* use the directly determined last velocity, not actually the averaged half steps */
1726 if (bTrotter && ir->eI == eiVV)
1728 enerd->term[F_EKIN] = last_ekin;
1730 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1734 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1738 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1740 /* ######### END PREPARING EDR OUTPUT ########### */
1745 gmx_bool do_dr, do_or;
1747 if (fplog && do_log && bDoExpanded)
1749 /* only needed if doing expanded ensemble */
1750 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1751 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1753 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1757 upd_mdebin(mdebin, bDoDHDL, TRUE,
1758 t, mdatoms->tmass, enerd, state,
1759 ir->fepvals, ir->expandedvals, lastbox,
1760 shake_vir, force_vir, total_vir, pres,
1761 ekind, mu_tot, constr);
1765 upd_mdebin_step(mdebin);
1768 do_dr = do_per_step(step, ir->nstdisreout);
1769 do_or = do_per_step(step, ir->nstorireout);
1771 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1773 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1775 if (ir->ePull != epullNO)
1777 pull_print_output(ir->pull, step, t);
1780 if (do_per_step(step, ir->nstlog))
1782 if (fflush(fplog) != 0)
1784 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1790 /* Have to do this part _after_ outputting the logfile and the edr file */
1791 /* Gets written into the state at the beginning of next loop*/
1792 state->fep_state = lamnew;
1794 /* Print the remaining wall clock time for the run */
1795 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1799 fprintf(stderr, "\n");
1801 print_time(stderr, walltime_accounting, step, ir, cr);
1804 /* Ion/water position swapping.
1805 * Not done in last step since trajectory writing happens before this call
1806 * in the MD loop and exchanges would be lost anyway. */
1807 bNeedRepartition = FALSE;
1808 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1809 do_per_step(step, ir->swap->nstswap))
1811 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1812 bRerunMD ? rerun_fr.x : state->x,
1813 bRerunMD ? rerun_fr.box : state->box,
1814 top_global, MASTER(cr) && bVerbose, bRerunMD);
1816 if (bNeedRepartition && DOMAINDECOMP(cr))
1818 dd_collect_state(cr->dd, state, state_global);
1822 /* Replica exchange */
1826 bExchanged = replica_exchange(fplog, cr, repl_ex,
1827 state_global, enerd,
1831 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1833 dd_partition_system(fplog, step, cr, TRUE, 1,
1834 state_global, top_global, ir,
1835 state, &f, mdatoms, top, fr,
1836 vsite, shellfc, constr,
1837 nrnb, wcycle, FALSE);
1842 bStartingFromCpt = FALSE;
1844 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1845 /* With all integrators, except VV, we need to retain the pressure
1846 * at the current step for coupling at the next step.
1848 if ((state->flags & (1<<estPRES_PREV)) &&
1850 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1852 /* Store the pressure in t_state for pressure coupling
1853 * at the next MD step.
1855 copy_mat(pres, state->pres_prev);
1858 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1860 if ( (membed != NULL) && (!bLastStep) )
1862 rescale_membed(step_rel, membed, state_global->x);
1869 /* read next frame from input trajectory */
1870 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1875 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1879 if (!bRerunMD || !rerun_fr.bStep)
1881 /* increase the MD step number */
1886 cycles = wallcycle_stop(wcycle, ewcSTEP);
1887 if (DOMAINDECOMP(cr) && wcycle)
1889 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1892 if (bPMETuneRunning || bPMETuneTry)
1894 /* PME grid + cut-off optimization with GPUs or PME nodes */
1896 /* Count the total cycles over the last steps */
1897 cycles_pmes += cycles;
1899 /* We can only switch cut-off at NS steps */
1900 if (step % ir->nstlist == 0)
1902 /* PME grid + cut-off optimization with GPUs or PME nodes */
1905 if (DDMASTER(cr->dd))
1907 /* PME node load is too high, start tuning */
1908 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1910 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1912 if (bPMETuneRunning &&
1913 fr->nbv->bUseGPU && DOMAINDECOMP(cr) &&
1914 !(cr->duty & DUTY_PME))
1916 /* Lock DLB=auto to off (does nothing when DLB=yes/no).
1917 * With GPUs + separate PME ranks, we don't want DLB.
1918 * This could happen when we scan coarse grids and
1919 * it would then never be turned off again.
1920 * This would hurt performance at the final, optimal
1921 * grid spacing, where DLB almost never helps.
1922 * Also, DLB can limit the cut-off for PME tuning.
1924 dd_dlb_set_lock(cr->dd, TRUE);
1927 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1929 bPMETuneTry = FALSE;
1932 if (bPMETuneRunning)
1934 /* init_step might not be a multiple of nstlist,
1935 * but the first cycle is always skipped anyhow.
1938 pme_load_balance(pme_loadbal, cr,
1939 (bVerbose && MASTER(cr)) ? stderr : NULL,
1941 ir, state, cycles_pmes,
1942 fr->ic, fr->nbv, &fr->pmedata,
1945 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1946 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1947 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1948 fr->rlist = fr->ic->rlist;
1949 fr->rlistlong = fr->ic->rlistlong;
1950 fr->rcoulomb = fr->ic->rcoulomb;
1951 fr->rvdw = fr->ic->rvdw;
1953 if (ir->eDispCorr != edispcNO)
1955 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1958 if (!bPMETuneRunning &&
1960 dd_dlb_is_locked(cr->dd))
1962 /* Unlock the DLB=auto, DLB is allowed to activate
1963 * (but we don't expect it to activate in most cases).
1965 dd_dlb_set_lock(cr->dd, FALSE);
1972 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1973 gs.set[eglsRESETCOUNTERS] != 0)
1975 /* Reset all the counters related to performance over the run */
1976 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1977 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1978 wcycle_set_reset_counters(wcycle, -1);
1979 if (!(cr->duty & DUTY_PME))
1981 /* Tell our PME node to reset its counters */
1982 gmx_pme_send_resetcounters(cr, step);
1984 /* Correct max_hours for the elapsed time */
1985 max_hours -= elapsed_time/(60.0*60.0);
1986 bResetCountersHalfMaxH = FALSE;
1987 gs.set[eglsRESETCOUNTERS] = 0;
1990 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1991 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1994 /* End of main MD loop */
1997 /* Closing TNG files can include compressing data. Therefore it is good to do that
1998 * before stopping the time measurements. */
1999 mdoutf_tng_close(outf);
2001 /* Stop measuring walltime */
2002 walltime_accounting_end(walltime_accounting);
2004 if (bRerunMD && MASTER(cr))
2009 if (!(cr->duty & DUTY_PME))
2011 /* Tell the PME only node to finish */
2012 gmx_pme_send_finish(cr);
2017 if (ir->nstcalcenergy > 0 && !bRerunMD)
2019 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
2020 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
2027 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2029 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)));
2030 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
2033 if (pme_loadbal != NULL)
2035 pme_loadbal_done(pme_loadbal, cr, fplog,
2036 fr->nbv != NULL && fr->nbv->bUseGPU);
2039 if (shellfc && fplog)
2041 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
2042 (nconverged*100.0)/step_rel);
2043 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
2047 if (repl_ex_nst > 0 && MASTER(cr))
2049 print_replica_exchange_statistics(fplog, repl_ex);
2052 /* IMD cleanup, if bIMD is TRUE. */
2053 IMD_finalize(ir->bIMD, ir->imd);
2055 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);