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, 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);
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, fr->cutoff_scheme == ecutsVERLET,
332 top_global, n_flexible_constraints(constr),
333 (ir->bContinuation ||
334 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
335 NULL : state_global->x);
337 if (shellfc && ir->eI == eiNM)
339 /* Currently shells don't work with Normal Modes */
340 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");
343 if (vsite && ir->eI == eiNM)
345 /* Currently virtual sites don't work with Normal Modes */
346 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");
351 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
352 set_deform_reference_box(upd,
353 deform_init_init_step_tpx,
354 deform_init_box_tpx);
355 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
359 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
360 if ((io > 2000) && MASTER(cr))
363 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
368 if (DOMAINDECOMP(cr))
370 top = dd_init_local_top(top_global);
373 dd_init_local_state(cr->dd, state_global, state);
375 if (DDMASTER(cr->dd) && ir->nstfout)
377 snew(f_global, state_global->natoms);
382 top = gmx_mtop_generate_local_top(top_global, ir);
384 forcerec_set_excl_load(fr, top);
386 state = serial_init_local_state(state_global);
389 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
393 set_vsite_top(vsite, top, mdatoms, cr);
396 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
398 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
403 make_local_shells(cr, mdatoms, shellfc);
406 setup_bonded_threading(fr, &top->idef);
409 /* Set up interactive MD (IMD) */
410 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
411 nfile, fnm, oenv, imdport, Flags);
413 if (DOMAINDECOMP(cr))
415 /* Distribute the charge groups over the nodes from the master node */
416 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
417 state_global, top_global, ir,
418 state, &f, mdatoms, top, fr,
419 vsite, shellfc, constr,
420 nrnb, wcycle, FALSE);
424 update_mdatoms(mdatoms, state->lambda[efptMASS]);
426 if (opt2bSet("-cpi", nfile, fnm))
428 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
432 bStateFromCP = FALSE;
437 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
444 /* Update mdebin with energy history if appending to output files */
445 if (Flags & MD_APPENDFILES)
447 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
451 /* We might have read an energy history from checkpoint,
452 * free the allocated memory and reset the counts.
454 done_energyhistory(&state_global->enerhist);
455 init_energyhistory(&state_global->enerhist);
458 /* Set the initial energy history in state by updating once */
459 update_energyhistory(&state_global->enerhist, mdebin);
462 /* Initialize constraints */
463 if (constr && !DOMAINDECOMP(cr))
465 set_constraints(constr, top, ir, mdatoms, cr);
470 /* We need to be sure replica exchange can only occur
471 * when the energies are current */
472 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
473 "repl_ex_nst", &repl_ex_nst);
474 /* This check needs to happen before inter-simulation
475 * signals are initialized, too */
477 if (repl_ex_nst > 0 && MASTER(cr))
479 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
480 repl_ex_nst, repl_ex_nex, repl_ex_seed);
483 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
484 * PME tuning is not supported with PME only for LJ and not for Coulomb.
486 if ((Flags & MD_TUNEPME) &&
487 EEL_PME(fr->eeltype) &&
488 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
491 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
493 if (cr->duty & DUTY_PME)
495 /* Start tuning right away, as we can't measure the load */
496 bPMETuneRunning = TRUE;
500 /* Separate PME nodes, we can measure the PP/PME load balance */
505 if (!ir->bContinuation && !bRerunMD)
507 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
509 /* Set the velocities of frozen particles to zero */
510 for (i = 0; i < mdatoms->homenr; i++)
512 for (m = 0; m < DIM; m++)
514 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
524 /* Constrain the initial coordinates and velocities */
525 do_constrain_first(fplog, constr, ir, mdatoms, state,
530 /* Construct the virtual sites for the initial configuration */
531 construct_vsites(vsite, state->x, ir->delta_t, NULL,
532 top->idef.iparams, top->idef.il,
533 fr->ePBC, fr->bMolPBC, cr, state->box);
539 /* set free energy calculation frequency as the minimum
540 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
541 nstfep = ir->fepvals->nstdhdl;
544 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
548 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
551 /* I'm assuming we need global communication the first time! MRS */
552 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
553 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
554 | (bVV ? CGLO_PRESSURE : 0)
555 | (bVV ? CGLO_CONSTRAINT : 0)
556 | (bRerunMD ? CGLO_RERUNMD : 0)
557 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
559 bSumEkinhOld = FALSE;
560 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
561 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
562 constr, NULL, FALSE, state->box,
563 top_global, &bSumEkinhOld, cglo_flags);
564 if (ir->eI == eiVVAK)
566 /* a second call to get the half step temperature initialized as well */
567 /* we do the same call as above, but turn the pressure off -- internally to
568 compute_globals, this is recognized as a velocity verlet half-step
569 kinetic energy calculation. This minimized excess variables, but
570 perhaps loses some logic?*/
572 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
573 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
574 constr, NULL, FALSE, state->box,
575 top_global, &bSumEkinhOld,
576 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
579 /* Calculate the initial half step temperature, and save the ekinh_old */
580 if (!(Flags & MD_STARTFROMCPT))
582 for (i = 0; (i < ir->opts.ngtc); i++)
584 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
589 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
590 and there is no previous step */
593 /* if using an iterative algorithm, we need to create a working directory for the state. */
596 bufstate = init_bufstate(state);
599 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
600 temperature control */
601 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
605 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
608 "RMS relative constraint deviation after constraining: %.2e\n",
609 constr_rmsd(constr, FALSE));
611 if (EI_STATE_VELOCITY(ir->eI))
613 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
617 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
618 " input trajectory '%s'\n\n",
619 *(top_global->name), opt2fn("-rerun", nfile, fnm));
622 fprintf(stderr, "Calculated time to finish depends on nsteps from "
623 "run input file,\nwhich may not correspond to the time "
624 "needed to process input trajectory.\n\n");
630 fprintf(stderr, "starting mdrun '%s'\n",
631 *(top_global->name));
634 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
638 sprintf(tbuf, "%s", "infinite");
640 if (ir->init_step > 0)
642 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
643 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
644 gmx_step_str(ir->init_step, sbuf2),
645 ir->init_step*ir->delta_t);
649 fprintf(stderr, "%s steps, %s ps.\n",
650 gmx_step_str(ir->nsteps, sbuf), tbuf);
653 fprintf(fplog, "\n");
656 walltime_accounting_start(walltime_accounting);
657 wallcycle_start(wcycle, ewcRUN);
658 print_start(fplog, cr, walltime_accounting, "mdrun");
660 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
662 chkpt_ret = fcCheckPointParallel( cr->nodeid,
666 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
671 /***********************************************************
675 ************************************************************/
677 /* if rerunMD then read coordinates and velocities from input trajectory */
680 if (getenv("GMX_FORCE_UPDATE"))
688 bNotLastFrame = read_first_frame(oenv, &status,
689 opt2fn("-rerun", nfile, fnm),
690 &rerun_fr, TRX_NEED_X | TRX_READ_V);
691 if (rerun_fr.natoms != top_global->natoms)
694 "Number of atoms in trajectory (%d) does not match the "
695 "run input file (%d)\n",
696 rerun_fr.natoms, top_global->natoms);
698 if (ir->ePBC != epbcNONE)
702 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);
704 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
706 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
713 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
716 if (ir->ePBC != epbcNONE)
718 /* Set the shift vectors.
719 * Necessary here when have a static box different from the tpr box.
721 calc_shifts(rerun_fr.box, fr->shift_vec);
725 /* loop over MD steps or if rerunMD to end of input trajectory */
727 /* Skip the first Nose-Hoover integration when we get the state from tpx */
728 bStateFromTPX = !bStateFromCP;
729 bInitStep = bFirstStep && (bStateFromTPX || bVV);
730 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
732 bSumEkinhOld = FALSE;
734 bNeedRepartition = FALSE;
736 init_global_signals(&gs, cr, ir, repl_ex_nst);
738 step = ir->init_step;
741 if (ir->nstlist == -1)
743 init_nlistheuristics(&nlh, bGStatEveryStep, step);
746 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
748 /* check how many steps are left in other sims */
749 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
753 /* and stop now if we should */
754 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
755 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
756 while (!bLastStep || (bRerunMD && bNotLastFrame))
759 wallcycle_start(wcycle, ewcSTEP);
765 step = rerun_fr.step;
766 step_rel = step - ir->init_step;
779 bLastStep = (step_rel == ir->nsteps);
780 t = t0 + step*ir->delta_t;
783 if (ir->efep != efepNO || ir->bSimTemp)
785 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
786 requiring different logic. */
788 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
789 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
790 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
791 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
792 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
797 update_annealing_target_temp(&(ir->opts), t);
802 if (!DOMAINDECOMP(cr) || MASTER(cr))
804 for (i = 0; i < state_global->natoms; i++)
806 copy_rvec(rerun_fr.x[i], state_global->x[i]);
810 for (i = 0; i < state_global->natoms; i++)
812 copy_rvec(rerun_fr.v[i], state_global->v[i]);
817 for (i = 0; i < state_global->natoms; i++)
819 clear_rvec(state_global->v[i]);
823 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
824 " Ekin, temperature and pressure are incorrect,\n"
825 " the virial will be incorrect when constraints are present.\n"
827 bRerunWarnNoV = FALSE;
831 copy_mat(rerun_fr.box, state_global->box);
832 copy_mat(state_global->box, state->box);
834 if (vsite && (Flags & MD_RERUN_VSITE))
836 if (DOMAINDECOMP(cr))
838 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
842 /* Following is necessary because the graph may get out of sync
843 * with the coordinates if we only have every N'th coordinate set
845 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
846 shift_self(graph, state->box, state->x);
848 construct_vsites(vsite, state->x, ir->delta_t, state->v,
849 top->idef.iparams, top->idef.il,
850 fr->ePBC, fr->bMolPBC, cr, state->box);
853 unshift_self(graph, state->box, state->x);
858 /* Stop Center of Mass motion */
859 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
863 /* for rerun MD always do Neighbour Searching */
864 bNS = (bFirstStep || ir->nstlist != 0);
869 /* Determine whether or not to do Neighbour Searching and LR */
870 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
872 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
873 (ir->nstlist == -1 && nlh.nabnsb > 0));
875 if (bNS && ir->nstlist == -1)
877 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
881 /* check whether we should stop because another simulation has
885 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
886 (multisim_nsteps != ir->nsteps) )
893 "Stopping simulation %d because another one has finished\n",
897 gs.sig[eglsCHKPT] = 1;
902 /* < 0 means stop at next step, > 0 means stop at next NS step */
903 if ( (gs.set[eglsSTOPCOND] < 0) ||
904 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
909 /* Determine whether or not to update the Born radii if doing GB */
910 bBornRadii = bFirstStep;
911 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
916 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
917 do_verbose = bVerbose &&
918 (step % stepout == 0 || bFirstStep || bLastStep);
920 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
928 bMasterState = FALSE;
929 /* Correct the new box if it is too skewed */
930 if (DYNAMIC_BOX(*ir))
932 if (correct_box(fplog, step, state->box, graph))
937 if (DOMAINDECOMP(cr) && bMasterState)
939 dd_collect_state(cr->dd, state, state_global);
943 if (DOMAINDECOMP(cr))
945 /* Repartition the domain decomposition */
946 wallcycle_start(wcycle, ewcDOMDEC);
947 dd_partition_system(fplog, step, cr,
948 bMasterState, nstglobalcomm,
949 state_global, top_global, ir,
950 state, &f, mdatoms, top, fr,
951 vsite, shellfc, constr,
953 do_verbose && !bPMETuneRunning);
954 wallcycle_stop(wcycle, ewcDOMDEC);
955 /* If using an iterative integrator, reallocate space to match the decomposition */
959 if (MASTER(cr) && do_log)
961 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
964 if (ir->efep != efepNO)
966 update_mdatoms(mdatoms, state->lambda[efptMASS]);
969 if ((bRerunMD && rerun_fr.bV) || bExchanged)
972 /* We need the kinetic energy at minus the half step for determining
973 * the full step kinetic energy and possibly for T-coupling.*/
974 /* This may not be quite working correctly yet . . . . */
975 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
976 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
977 constr, NULL, FALSE, state->box,
978 top_global, &bSumEkinhOld,
979 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
981 clear_mat(force_vir);
983 /* We write a checkpoint at this MD step when:
984 * either at an NS step when we signalled through gs,
985 * or at the last step (but not when we do not want confout),
986 * but never at the first step or with rerun.
988 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
989 (bLastStep && (Flags & MD_CONFOUT))) &&
990 step > ir->init_step && !bRerunMD);
993 gs.set[eglsCHKPT] = 0;
996 /* Determine the energy and pressure:
997 * at nstcalcenergy steps and at energy output steps (set below).
999 if (EI_VV(ir->eI) && (!bInitStep))
1001 /* for vv, the first half of the integration actually corresponds
1002 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1003 but the virial needs to be calculated on both the current step and the 'next' step. Future
1004 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1006 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1007 bCalcVir = bCalcEner ||
1008 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1012 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1013 bCalcVir = bCalcEner ||
1014 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1017 /* Do we need global communication ? */
1018 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1019 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1020 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1022 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1024 if (do_ene || do_log)
1031 /* these CGLO_ options remain the same throughout the iteration */
1032 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1033 (bGStat ? CGLO_GSTAT : 0)
1036 force_flags = (GMX_FORCE_STATECHANGED |
1037 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1038 GMX_FORCE_ALLFORCES |
1040 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1041 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1042 (bDoFEP ? GMX_FORCE_DHDL : 0)
1047 if (do_per_step(step, ir->nstcalclr))
1049 force_flags |= GMX_FORCE_DO_LR;
1055 /* Now is the time to relax the shells */
1056 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1057 ir, bNS, force_flags,
1060 state, f, force_vir, mdatoms,
1061 nrnb, wcycle, graph, groups,
1062 shellfc, fr, bBornRadii, t, mu_tot,
1064 mdoutf_get_fp_field(outf));
1074 /* The coordinates (x) are shifted (to get whole molecules)
1076 * This is parallellized as well, and does communication too.
1077 * Check comments in sim_util.c
1079 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1080 state->box, state->x, &state->hist,
1081 f, force_vir, mdatoms, enerd, fcd,
1082 state->lambda, graph,
1083 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1084 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1087 if (bVV && !bStartingFromCpt && !bRerunMD)
1088 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1090 if (ir->eI == eiVV && bInitStep)
1092 /* if using velocity verlet with full time step Ekin,
1093 * take the first half step only to compute the
1094 * virial for the first step. From there,
1095 * revert back to the initial coordinates
1096 * so that the input is actually the initial step.
1098 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1102 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1103 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1106 /* If we are using twin-range interactions where the long-range component
1107 * is only evaluated every nstcalclr>1 steps, we should do a special update
1108 * step to combine the long-range forces on these steps.
1109 * For nstcalclr=1 this is not done, since the forces would have been added
1110 * directly to the short-range forces already.
1112 * TODO Remove various aspects of VV+twin-range in master
1113 * branch, because VV integrators did not ever support
1114 * twin-range multiple time stepping with constraints.
1116 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1118 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1119 f, bUpdateDoLR, fr->f_twin, fcd,
1120 ekind, M, upd, bInitStep, etrtVELOCITY1,
1121 cr, nrnb, constr, &top->idef);
1123 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1125 gmx_iterate_init(&iterate, TRUE);
1127 /* for iterations, we save these vectors, as we will be self-consistently iterating
1130 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1132 /* save the state */
1133 if (iterate.bIterationActive)
1135 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1138 bFirstIterate = TRUE;
1139 while (bFirstIterate || iterate.bIterationActive)
1141 if (iterate.bIterationActive)
1143 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1144 if (bFirstIterate && bTrotter)
1146 /* The first time through, we need a decent first estimate
1147 of veta(t+dt) to compute the constraints. Do
1148 this by computing the box volume part of the
1149 trotter integration at this time. Nothing else
1150 should be changed by this routine here. If
1151 !(first time), we start with the previous value
1154 veta_save = state->veta;
1155 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1156 vetanew = state->veta;
1157 state->veta = veta_save;
1162 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1164 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1165 state, fr->bMolPBC, graph, f,
1166 &top->idef, shake_vir,
1167 cr, nrnb, wcycle, upd, constr,
1168 TRUE, bCalcVir, vetanew);
1172 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1178 /* Need to unshift here if a do_force has been
1179 called in the previous step */
1180 unshift_self(graph, state->box, state->x);
1183 /* if VV, compute the pressure and constraints */
1184 /* For VV2, we strictly only need this if using pressure
1185 * control, but we really would like to have accurate pressures
1187 * Think about ways around this in the future?
1188 * For now, keep this choice in comments.
1190 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1191 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1193 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1194 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1196 bSumEkinhOld = TRUE;
1198 /* for vv, the first half of the integration actually corresponds to the previous step.
1199 So we need information from the last step in the first half of the integration */
1200 if (bGStat || do_per_step(step-1, nstglobalcomm))
1202 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1203 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1204 constr, NULL, FALSE, state->box,
1205 top_global, &bSumEkinhOld,
1208 | (bTemp ? CGLO_TEMPERATURE : 0)
1209 | (bPres ? CGLO_PRESSURE : 0)
1210 | (bPres ? CGLO_CONSTRAINT : 0)
1211 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1212 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1215 /* explanation of above:
1216 a) We compute Ekin at the full time step
1217 if 1) we are using the AveVel Ekin, and it's not the
1218 initial step, or 2) if we are using AveEkin, but need the full
1219 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1220 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1221 EkinAveVel because it's needed for the pressure */
1223 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1228 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1229 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1236 /* We need the kinetic energy at minus the half step for determining
1237 * the full step kinetic energy and possibly for T-coupling.*/
1238 /* This may not be quite working correctly yet . . . . */
1239 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1240 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1241 constr, NULL, FALSE, state->box,
1242 top_global, &bSumEkinhOld,
1243 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1248 if (iterate.bIterationActive &&
1249 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1250 state->veta, &vetanew))
1254 bFirstIterate = FALSE;
1257 if (bTrotter && !bInitStep)
1259 copy_mat(shake_vir, state->svir_prev);
1260 copy_mat(force_vir, state->fvir_prev);
1261 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1263 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1264 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1265 enerd->term[F_EKIN] = trace(ekind->ekin);
1268 /* if it's the initial step, we performed this first step just to get the constraint virial */
1269 if (bInitStep && ir->eI == eiVV)
1271 copy_rvecn(cbuf, state->v, 0, state->natoms);
1275 /* MRS -- now done iterating -- compute the conserved quantity */
1278 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1281 last_ekin = enerd->term[F_EKIN];
1283 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1285 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1287 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1290 sum_dhdl(enerd, state->lambda, ir->fepvals);
1294 /* ######## END FIRST UPDATE STEP ############## */
1295 /* ######## If doing VV, we now have v(dt) ###### */
1298 /* perform extended ensemble sampling in lambda - we don't
1299 actually move to the new state before outputting
1300 statistics, but if performing simulated tempering, we
1301 do update the velocities and the tau_t. */
1303 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1304 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1305 copy_df_history(&state_global->dfhist, &state->dfhist);
1308 /* Now we have the energies and forces corresponding to the
1309 * coordinates at time t. We must output all of this before
1312 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1313 ir, state, state_global, top_global, fr,
1314 outf, mdebin, ekind, f, f_global,
1316 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1318 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1319 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1321 /* kludge -- virial is lost with restart for NPT control. Must restart */
1322 if (bStartingFromCpt && bVV)
1324 copy_mat(state->svir_prev, shake_vir);
1325 copy_mat(state->fvir_prev, force_vir);
1328 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1330 /* Check whether everything is still allright */
1331 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1332 #ifdef GMX_THREAD_MPI
1337 /* this is just make gs.sig compatible with the hack
1338 of sending signals around by MPI_Reduce with together with
1340 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1342 gs.sig[eglsSTOPCOND] = 1;
1344 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1346 gs.sig[eglsSTOPCOND] = -1;
1348 /* < 0 means stop at next step, > 0 means stop at next NS step */
1352 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1353 gmx_get_signal_name(),
1354 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1358 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1359 gmx_get_signal_name(),
1360 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1362 handled_stop_condition = (int)gmx_get_stop_condition();
1364 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1365 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1366 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1368 /* Signal to terminate the run */
1369 gs.sig[eglsSTOPCOND] = 1;
1372 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1374 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1377 if (bResetCountersHalfMaxH && MASTER(cr) &&
1378 elapsed_time > max_hours*60.0*60.0*0.495)
1380 gs.sig[eglsRESETCOUNTERS] = 1;
1383 if (ir->nstlist == -1 && !bRerunMD)
1385 /* When bGStatEveryStep=FALSE, global_stat is only called
1386 * when we check the atom displacements, not at NS steps.
1387 * This means that also the bonded interaction count check is not
1388 * performed immediately after NS. Therefore a few MD steps could
1389 * be performed with missing interactions.
1390 * But wrong energies are never written to file,
1391 * since energies are only written after global_stat
1394 if (step >= nlh.step_nscheck)
1396 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1397 nlh.scale_tot, state->x);
1401 /* This is not necessarily true,
1402 * but step_nscheck is determined quite conservatively.
1408 /* In parallel we only have to check for checkpointing in steps
1409 * where we do global communication,
1410 * otherwise the other nodes don't know.
1412 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1415 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1416 gs.set[eglsCHKPT] == 0)
1418 gs.sig[eglsCHKPT] = 1;
1421 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1426 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1428 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1430 gmx_bool bIfRandomize;
1431 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1432 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1433 if (constr && bIfRandomize)
1435 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1436 state, fr->bMolPBC, graph, f,
1437 &top->idef, tmp_vir,
1438 cr, nrnb, wcycle, upd, constr,
1439 TRUE, bCalcVir, vetanew);
1444 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1446 gmx_iterate_init(&iterate, TRUE);
1447 /* for iterations, we save these vectors, as we will be redoing the calculations */
1448 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1451 bFirstIterate = TRUE;
1452 while (bFirstIterate || iterate.bIterationActive)
1454 /* We now restore these vectors to redo the calculation with improved extended variables */
1455 if (iterate.bIterationActive)
1457 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1460 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1461 so scroll down for that logic */
1463 /* ######### START SECOND UPDATE STEP ################# */
1464 /* Box is changed in update() when we do pressure coupling,
1465 * but we should still use the old box for energy corrections and when
1466 * writing it to the energy file, so it matches the trajectory files for
1467 * the same timestep above. Make a copy in a separate array.
1469 copy_mat(state->box, lastbox);
1474 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1476 wallcycle_start(wcycle, ewcUPDATE);
1477 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1480 if (iterate.bIterationActive)
1488 /* we use a new value of scalevir to converge the iterations faster */
1489 scalevir = tracevir/trace(shake_vir);
1491 msmul(shake_vir, scalevir, shake_vir);
1492 m_add(force_vir, shake_vir, total_vir);
1493 clear_mat(shake_vir);
1495 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1496 /* We can only do Berendsen coupling after we have summed
1497 * the kinetic energy or virial. Since the happens
1498 * in global_state after update, we should only do it at
1499 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1504 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1505 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1510 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1512 /* velocity half-step update */
1513 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1514 bUpdateDoLR, fr->f_twin, fcd,
1515 ekind, M, upd, FALSE, etrtVELOCITY2,
1516 cr, nrnb, constr, &top->idef);
1519 /* Above, initialize just copies ekinh into ekin,
1520 * it doesn't copy position (for VV),
1521 * and entire integrator for MD.
1524 if (ir->eI == eiVVAK)
1526 copy_rvecn(state->x, cbuf, 0, state->natoms);
1528 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1530 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1531 bUpdateDoLR, fr->f_twin, fcd,
1532 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1533 wallcycle_stop(wcycle, ewcUPDATE);
1535 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1536 fr->bMolPBC, graph, f,
1537 &top->idef, shake_vir,
1538 cr, nrnb, wcycle, upd, constr,
1539 FALSE, bCalcVir, state->veta);
1541 if (ir->eI == eiVVAK)
1543 /* erase F_EKIN and F_TEMP here? */
1544 /* just compute the kinetic energy at the half step to perform a trotter step */
1545 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1546 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1547 constr, NULL, FALSE, lastbox,
1548 top_global, &bSumEkinhOld,
1549 cglo_flags | CGLO_TEMPERATURE
1551 wallcycle_start(wcycle, ewcUPDATE);
1552 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1553 /* now we know the scaling, we can compute the positions again again */
1554 copy_rvecn(cbuf, state->x, 0, state->natoms);
1556 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1558 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1559 bUpdateDoLR, fr->f_twin, fcd,
1560 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1561 wallcycle_stop(wcycle, ewcUPDATE);
1563 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1564 /* are the small terms in the shake_vir here due
1565 * to numerical errors, or are they important
1566 * physically? I'm thinking they are just errors, but not completely sure.
1567 * For now, will call without actually constraining, constr=NULL*/
1568 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1569 state, fr->bMolPBC, graph, f,
1570 &top->idef, tmp_vir,
1571 cr, nrnb, wcycle, upd, NULL,
1577 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1580 if (fr->bSepDVDL && fplog && do_log)
1582 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1586 /* this factor or 2 correction is necessary
1587 because half of the constraint force is removed
1588 in the vv step, so we have to double it. See
1589 the Redmine issue #1255. It is not yet clear
1590 if the factor of 2 is exact, or just a very
1591 good approximation, and this will be
1592 investigated. The next step is to see if this
1593 can be done adding a dhdl contribution from the
1594 rattle step, but this is somewhat more
1595 complicated with the current code. Will be
1596 investigated, hopefully for 4.6.3. However,
1597 this current solution is much better than
1598 having it completely wrong.
1600 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1604 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1609 /* Need to unshift here */
1610 unshift_self(graph, state->box, state->x);
1615 wallcycle_start(wcycle, ewcVSITECONSTR);
1618 shift_self(graph, state->box, state->x);
1620 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1621 top->idef.iparams, top->idef.il,
1622 fr->ePBC, fr->bMolPBC, cr, state->box);
1626 unshift_self(graph, state->box, state->x);
1628 wallcycle_stop(wcycle, ewcVSITECONSTR);
1631 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1632 /* With Leap-Frog we can skip compute_globals at
1633 * non-communication steps, but we need to calculate
1634 * the kinetic energy one step before communication.
1636 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1638 if (ir->nstlist == -1 && bFirstIterate)
1640 gs.sig[eglsNABNSB] = nlh.nabnsb;
1642 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1643 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1645 bFirstIterate ? &gs : NULL,
1646 (step_rel % gs.nstms == 0) &&
1647 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1649 top_global, &bSumEkinhOld,
1651 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1652 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1653 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1654 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1655 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1656 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1659 if (ir->nstlist == -1 && bFirstIterate)
1661 nlh.nabnsb = gs.set[eglsNABNSB];
1662 gs.set[eglsNABNSB] = 0;
1665 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1666 /* ############# END CALC EKIN AND PRESSURE ################# */
1668 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1669 the virial that should probably be addressed eventually. state->veta has better properies,
1670 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1671 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1673 if (iterate.bIterationActive &&
1674 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1675 trace(shake_vir), &tracevir))
1679 bFirstIterate = FALSE;
1682 if (!bVV || bRerunMD)
1684 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1685 sum_dhdl(enerd, state->lambda, ir->fepvals);
1687 update_box(fplog, step, ir, mdatoms, state, f,
1688 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1690 /* ################# END UPDATE STEP 2 ################# */
1691 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1693 /* The coordinates (x) were unshifted in update */
1696 /* We will not sum ekinh_old,
1697 * so signal that we still have to do it.
1699 bSumEkinhOld = TRUE;
1702 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1704 /* use the directly determined last velocity, not actually the averaged half steps */
1705 if (bTrotter && ir->eI == eiVV)
1707 enerd->term[F_EKIN] = last_ekin;
1709 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1713 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1717 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1719 /* ######### END PREPARING EDR OUTPUT ########### */
1724 gmx_bool do_dr, do_or;
1726 if (fplog && do_log && bDoExpanded)
1728 /* only needed if doing expanded ensemble */
1729 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1730 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1732 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1736 upd_mdebin(mdebin, bDoDHDL, TRUE,
1737 t, mdatoms->tmass, enerd, state,
1738 ir->fepvals, ir->expandedvals, lastbox,
1739 shake_vir, force_vir, total_vir, pres,
1740 ekind, mu_tot, constr);
1744 upd_mdebin_step(mdebin);
1747 do_dr = do_per_step(step, ir->nstdisreout);
1748 do_or = do_per_step(step, ir->nstorireout);
1750 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1752 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1754 if (ir->ePull != epullNO)
1756 pull_print_output(ir->pull, step, t);
1759 if (do_per_step(step, ir->nstlog))
1761 if (fflush(fplog) != 0)
1763 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1769 /* Have to do this part _after_ outputting the logfile and the edr file */
1770 /* Gets written into the state at the beginning of next loop*/
1771 state->fep_state = lamnew;
1773 /* Print the remaining wall clock time for the run */
1774 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1778 fprintf(stderr, "\n");
1780 print_time(stderr, walltime_accounting, step, ir, cr);
1783 /* Ion/water position swapping.
1784 * Not done in last step since trajectory writing happens before this call
1785 * in the MD loop and exchanges would be lost anyway. */
1786 bNeedRepartition = FALSE;
1787 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1788 do_per_step(step, ir->swap->nstswap))
1790 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1791 bRerunMD ? rerun_fr.x : state->x,
1792 bRerunMD ? rerun_fr.box : state->box,
1793 top_global, MASTER(cr) && bVerbose, bRerunMD);
1795 if (bNeedRepartition && DOMAINDECOMP(cr))
1797 dd_collect_state(cr->dd, state, state_global);
1801 /* Replica exchange */
1803 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1804 do_per_step(step, repl_ex_nst))
1806 bExchanged = replica_exchange(fplog, cr, repl_ex,
1807 state_global, enerd,
1811 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1813 dd_partition_system(fplog, step, cr, TRUE, 1,
1814 state_global, top_global, ir,
1815 state, &f, mdatoms, top, fr,
1816 vsite, shellfc, constr,
1817 nrnb, wcycle, FALSE);
1822 bStartingFromCpt = FALSE;
1824 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1825 /* With all integrators, except VV, we need to retain the pressure
1826 * at the current step for coupling at the next step.
1828 if ((state->flags & (1<<estPRES_PREV)) &&
1830 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1832 /* Store the pressure in t_state for pressure coupling
1833 * at the next MD step.
1835 copy_mat(pres, state->pres_prev);
1838 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1840 if ( (membed != NULL) && (!bLastStep) )
1842 rescale_membed(step_rel, membed, state_global->x);
1849 /* read next frame from input trajectory */
1850 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1855 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1859 if (!bRerunMD || !rerun_fr.bStep)
1861 /* increase the MD step number */
1866 cycles = wallcycle_stop(wcycle, ewcSTEP);
1867 if (DOMAINDECOMP(cr) && wcycle)
1869 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1872 if (bPMETuneRunning || bPMETuneTry)
1874 /* PME grid + cut-off optimization with GPUs or PME nodes */
1876 /* Count the total cycles over the last steps */
1877 cycles_pmes += cycles;
1879 /* We can only switch cut-off at NS steps */
1880 if (step % ir->nstlist == 0)
1882 /* PME grid + cut-off optimization with GPUs or PME nodes */
1885 if (DDMASTER(cr->dd))
1887 /* PME node load is too high, start tuning */
1888 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1890 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1892 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1894 bPMETuneTry = FALSE;
1897 if (bPMETuneRunning)
1899 /* init_step might not be a multiple of nstlist,
1900 * but the first cycle is always skipped anyhow.
1903 pme_load_balance(pme_loadbal, cr,
1904 (bVerbose && MASTER(cr)) ? stderr : NULL,
1906 ir, state, cycles_pmes,
1907 fr->ic, fr->nbv, &fr->pmedata,
1910 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1911 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1912 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1913 fr->rlist = fr->ic->rlist;
1914 fr->rlistlong = fr->ic->rlistlong;
1915 fr->rcoulomb = fr->ic->rcoulomb;
1916 fr->rvdw = fr->ic->rvdw;
1918 if (ir->eDispCorr != edispcNO)
1920 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1927 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1928 gs.set[eglsRESETCOUNTERS] != 0)
1930 /* Reset all the counters related to performance over the run */
1931 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1932 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1933 wcycle_set_reset_counters(wcycle, -1);
1934 if (!(cr->duty & DUTY_PME))
1936 /* Tell our PME node to reset its counters */
1937 gmx_pme_send_resetcounters(cr, step);
1939 /* Correct max_hours for the elapsed time */
1940 max_hours -= elapsed_time/(60.0*60.0);
1941 bResetCountersHalfMaxH = FALSE;
1942 gs.set[eglsRESETCOUNTERS] = 0;
1945 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1946 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1949 /* End of main MD loop */
1952 /* Stop measuring walltime */
1953 walltime_accounting_end(walltime_accounting);
1955 if (bRerunMD && MASTER(cr))
1960 if (!(cr->duty & DUTY_PME))
1962 /* Tell the PME only node to finish */
1963 gmx_pme_send_finish(cr);
1968 if (ir->nstcalcenergy > 0 && !bRerunMD)
1970 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1971 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1978 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1980 fprintf(fplog, "Average neighborlist lifetime: %.1f steps, std.dev.: %.1f steps\n", nlh.s1/nlh.nns, sqrt(nlh.s2/nlh.nns - sqr(nlh.s1/nlh.nns)));
1981 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1984 if (pme_loadbal != NULL)
1986 pme_loadbal_done(pme_loadbal, cr, fplog,
1987 fr->nbv != NULL && fr->nbv->bUseGPU);
1990 if (shellfc && fplog)
1992 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1993 (nconverged*100.0)/step_rel);
1994 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1998 if (repl_ex_nst > 0 && MASTER(cr))
2000 print_replica_exchange_statistics(fplog, repl_ex);
2003 /* IMD cleanup, if bIMD is TRUE. */
2004 IMD_finalize(ir->bIMD, ir->imd);
2006 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);