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.
44 #include "gromacs/math/vec.h"
53 #include "md_support.h"
54 #include "md_logging.h"
66 #include "domdec_network.h"
70 #include "gromacs/gmxpreprocess/compute_io.h"
71 #include "checkpoint.h"
72 #include "gromacs/topology/mtop_util.h"
73 #include "sighandler.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "pme_loadbal.h"
79 #include "types/nlistheuristics.h"
80 #include "types/iteratedconstraints.h"
81 #include "nbnxn_cuda_data_mgmt.h"
83 #include "gromacs/fileio/confio.h"
84 #include "gromacs/fileio/mdoutf.h"
85 #include "gromacs/fileio/trajectory_writing.h"
86 #include "gromacs/fileio/trnio.h"
87 #include "gromacs/fileio/trxio.h"
88 #include "gromacs/fileio/xtcio.h"
89 #include "gromacs/imd/imd.h"
90 #include "gromacs/pbcutil/mshift.h"
91 #include "gromacs/pbcutil/pbc.h"
92 #include "gromacs/pulling/pull.h"
93 #include "gromacs/swap/swapcoords.h"
94 #include "gromacs/timing/wallcycle.h"
95 #include "gromacs/timing/walltime_accounting.h"
96 #include "gromacs/utility/gmxmpi.h"
97 #include "gromacs/utility/smalloc.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;
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;
193 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
194 gmx_bool bResetCountersHalfMaxH = FALSE;
195 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
196 gmx_bool bUpdateDoLR;
200 real veta_save, scalevir, tracevir;
206 real saved_conserved_quantity = 0;
210 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
211 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
212 gmx_iterate_t iterate;
213 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
214 simulation stops. If equal to zero, don't
215 communicate any more between multisims.*/
216 /* PME load balancing data for GPU kernels */
217 pme_load_balancing_t pme_loadbal = NULL;
219 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
222 gmx_bool bIMDstep = FALSE;
225 /* Temporary addition for FAHCORE checkpointing */
229 /* Check for special mdrun options */
230 bRerunMD = (Flags & MD_RERUN);
231 if (Flags & MD_RESETCOUNTERSHALFWAY)
235 /* Signal to reset the counters half the simulation steps. */
236 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
238 /* Signal to reset the counters halfway the simulation time. */
239 bResetCountersHalfMaxH = (max_hours > 0);
242 /* md-vv uses averaged full step velocities for T-control
243 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
244 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
246 if (bVV) /* to store the initial velocities while computing virial */
248 snew(cbuf, top_global->natoms);
250 /* all the iteratative cases - only if there are constraints */
251 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
252 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
253 false in this step. The correct value, true or false,
254 is set at each step, as it depends on the frequency of temperature
255 and pressure control.*/
256 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
260 /* Since we don't know if the frames read are related in any way,
261 * rebuild the neighborlist at every step.
264 ir->nstcalcenergy = 1;
268 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
270 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
271 bGStatEveryStep = (nstglobalcomm == 1);
273 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
276 "To reduce the energy communication with nstlist = -1\n"
277 "the neighbor list validity should not be checked at every step,\n"
278 "this means that exact integration is not guaranteed.\n"
279 "The neighbor list validity is checked after:\n"
280 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
281 "In most cases this will result in exact integration.\n"
282 "This reduces the energy communication by a factor of 2 to 3.\n"
283 "If you want less energy communication, set nstlist > 3.\n\n");
288 ir->nstxout_compressed = 0;
290 groups = &top_global->groups;
293 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
294 &(state_global->fep_state), lam0,
295 nrnb, top_global, &upd,
296 nfile, fnm, &outf, &mdebin,
297 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
299 clear_mat(total_vir);
301 /* Energy terms and groups */
303 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
305 if (DOMAINDECOMP(cr))
311 snew(f, top_global->natoms);
314 /* Kinetic energy data */
316 init_ekindata(fplog, top_global, &(ir->opts), ekind);
317 /* needed for iteration of constraints */
319 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
320 /* Copy the cos acceleration to the groups struct */
321 ekind->cosacc.cos_accel = ir->cos_accel;
323 gstat = global_stat_init(ir);
326 /* Check for polarizable models and flexible constraints */
327 shellfc = init_shell_flexcon(fplog, fr->cutoff_scheme == ecutsVERLET,
328 top_global, n_flexible_constraints(constr),
329 (ir->bContinuation ||
330 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
331 NULL : state_global->x);
333 if (shellfc && ir->eI == eiNM)
335 /* Currently shells don't work with Normal Modes */
336 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");
339 if (vsite && ir->eI == eiNM)
341 /* Currently virtual sites don't work with Normal Modes */
342 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");
347 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
348 set_deform_reference_box(upd,
349 deform_init_init_step_tpx,
350 deform_init_box_tpx);
351 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
355 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
356 if ((io > 2000) && MASTER(cr))
359 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
364 if (DOMAINDECOMP(cr))
366 top = dd_init_local_top(top_global);
369 dd_init_local_state(cr->dd, state_global, state);
371 if (DDMASTER(cr->dd) && ir->nstfout)
373 snew(f_global, state_global->natoms);
378 top = gmx_mtop_generate_local_top(top_global, ir);
380 forcerec_set_excl_load(fr, top);
382 state = serial_init_local_state(state_global);
385 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
389 set_vsite_top(vsite, top, mdatoms, cr);
392 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
394 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
399 make_local_shells(cr, mdatoms, shellfc);
402 setup_bonded_threading(fr, &top->idef);
405 /* Set up interactive MD (IMD) */
406 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
407 nfile, fnm, oenv, imdport, Flags);
409 if (DOMAINDECOMP(cr))
411 /* Distribute the charge groups over the nodes from the master node */
412 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
413 state_global, top_global, ir,
414 state, &f, mdatoms, top, fr,
415 vsite, shellfc, constr,
416 nrnb, wcycle, FALSE);
420 update_mdatoms(mdatoms, state->lambda[efptMASS]);
422 if (opt2bSet("-cpi", nfile, fnm))
424 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
428 bStateFromCP = FALSE;
433 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
440 /* Update mdebin with energy history if appending to output files */
441 if (Flags & MD_APPENDFILES)
443 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
447 /* We might have read an energy history from checkpoint,
448 * free the allocated memory and reset the counts.
450 done_energyhistory(&state_global->enerhist);
451 init_energyhistory(&state_global->enerhist);
454 /* Set the initial energy history in state by updating once */
455 update_energyhistory(&state_global->enerhist, mdebin);
458 /* Initialize constraints */
459 if (constr && !DOMAINDECOMP(cr))
461 set_constraints(constr, top, ir, mdatoms, cr);
464 if (repl_ex_nst > 0 && MASTER(cr))
466 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
467 repl_ex_nst, repl_ex_nex, repl_ex_seed);
470 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
471 * PME tuning is not supported with PME only for LJ and not for Coulomb.
473 if ((Flags & MD_TUNEPME) &&
474 EEL_PME(fr->eeltype) &&
475 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
478 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
480 if (cr->duty & DUTY_PME)
482 /* Start tuning right away, as we can't measure the load */
483 bPMETuneRunning = TRUE;
487 /* Separate PME nodes, we can measure the PP/PME load balance */
492 if (!ir->bContinuation && !bRerunMD)
494 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
496 /* Set the velocities of frozen particles to zero */
497 for (i = 0; i < mdatoms->homenr; i++)
499 for (m = 0; m < DIM; m++)
501 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
511 /* Constrain the initial coordinates and velocities */
512 do_constrain_first(fplog, constr, ir, mdatoms, state,
517 /* Construct the virtual sites for the initial configuration */
518 construct_vsites(vsite, state->x, ir->delta_t, NULL,
519 top->idef.iparams, top->idef.il,
520 fr->ePBC, fr->bMolPBC, cr, state->box);
526 /* set free energy calculation frequency as the minimum
527 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
528 nstfep = ir->fepvals->nstdhdl;
531 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
535 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
538 /* I'm assuming we need global communication the first time! MRS */
539 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
540 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
541 | (bVV ? CGLO_PRESSURE : 0)
542 | (bVV ? CGLO_CONSTRAINT : 0)
543 | (bRerunMD ? CGLO_RERUNMD : 0)
544 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
546 bSumEkinhOld = FALSE;
547 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
548 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
549 constr, NULL, FALSE, state->box,
550 top_global, &bSumEkinhOld, cglo_flags);
551 if (ir->eI == eiVVAK)
553 /* a second call to get the half step temperature initialized as well */
554 /* we do the same call as above, but turn the pressure off -- internally to
555 compute_globals, this is recognized as a velocity verlet half-step
556 kinetic energy calculation. This minimized excess variables, but
557 perhaps loses some logic?*/
559 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
560 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
561 constr, NULL, FALSE, state->box,
562 top_global, &bSumEkinhOld,
563 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
566 /* Calculate the initial half step temperature, and save the ekinh_old */
567 if (!(Flags & MD_STARTFROMCPT))
569 for (i = 0; (i < ir->opts.ngtc); i++)
571 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
576 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
577 and there is no previous step */
580 /* if using an iterative algorithm, we need to create a working directory for the state. */
583 bufstate = init_bufstate(state);
586 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
587 temperature control */
588 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
592 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
595 "RMS relative constraint deviation after constraining: %.2e\n",
596 constr_rmsd(constr, FALSE));
598 if (EI_STATE_VELOCITY(ir->eI))
600 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
604 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
605 " input trajectory '%s'\n\n",
606 *(top_global->name), opt2fn("-rerun", nfile, fnm));
609 fprintf(stderr, "Calculated time to finish depends on nsteps from "
610 "run input file,\nwhich may not correspond to the time "
611 "needed to process input trajectory.\n\n");
617 fprintf(stderr, "starting mdrun '%s'\n",
618 *(top_global->name));
621 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
625 sprintf(tbuf, "%s", "infinite");
627 if (ir->init_step > 0)
629 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
630 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
631 gmx_step_str(ir->init_step, sbuf2),
632 ir->init_step*ir->delta_t);
636 fprintf(stderr, "%s steps, %s ps.\n",
637 gmx_step_str(ir->nsteps, sbuf), tbuf);
640 fprintf(fplog, "\n");
643 walltime_accounting_start(walltime_accounting);
644 wallcycle_start(wcycle, ewcRUN);
645 print_start(fplog, cr, walltime_accounting, "mdrun");
647 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
649 chkpt_ret = fcCheckPointParallel( cr->nodeid,
653 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
658 /***********************************************************
662 ************************************************************/
664 /* if rerunMD then read coordinates and velocities from input trajectory */
667 if (getenv("GMX_FORCE_UPDATE"))
675 bNotLastFrame = read_first_frame(oenv, &status,
676 opt2fn("-rerun", nfile, fnm),
677 &rerun_fr, TRX_NEED_X | TRX_READ_V);
678 if (rerun_fr.natoms != top_global->natoms)
681 "Number of atoms in trajectory (%d) does not match the "
682 "run input file (%d)\n",
683 rerun_fr.natoms, top_global->natoms);
685 if (ir->ePBC != epbcNONE)
689 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);
691 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
693 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
700 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
703 if (ir->ePBC != epbcNONE)
705 /* Set the shift vectors.
706 * Necessary here when have a static box different from the tpr box.
708 calc_shifts(rerun_fr.box, fr->shift_vec);
712 /* loop over MD steps or if rerunMD to end of input trajectory */
714 /* Skip the first Nose-Hoover integration when we get the state from tpx */
715 bStateFromTPX = !bStateFromCP;
716 bInitStep = bFirstStep && (bStateFromTPX || bVV);
717 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
718 bSumEkinhOld = FALSE;
720 bNeedRepartition = FALSE;
722 init_global_signals(&gs, cr, ir, repl_ex_nst);
724 step = ir->init_step;
727 init_nlistheuristics(&nlh, bGStatEveryStep, step);
729 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
731 /* check how many steps are left in other sims */
732 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
736 /* and stop now if we should */
737 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
738 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
739 while (!bLastStep || (bRerunMD && bNotLastFrame))
742 wallcycle_start(wcycle, ewcSTEP);
748 step = rerun_fr.step;
749 step_rel = step - ir->init_step;
762 bLastStep = (step_rel == ir->nsteps);
763 t = t0 + step*ir->delta_t;
766 if (ir->efep != efepNO || ir->bSimTemp)
768 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
769 requiring different logic. */
771 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
772 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
773 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
774 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
775 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
778 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
779 do_per_step(step, repl_ex_nst));
783 update_annealing_target_temp(&(ir->opts), t);
788 if (!DOMAINDECOMP(cr) || MASTER(cr))
790 for (i = 0; i < state_global->natoms; i++)
792 copy_rvec(rerun_fr.x[i], state_global->x[i]);
796 for (i = 0; i < state_global->natoms; i++)
798 copy_rvec(rerun_fr.v[i], state_global->v[i]);
803 for (i = 0; i < state_global->natoms; i++)
805 clear_rvec(state_global->v[i]);
809 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
810 " Ekin, temperature and pressure are incorrect,\n"
811 " the virial will be incorrect when constraints are present.\n"
813 bRerunWarnNoV = FALSE;
817 copy_mat(rerun_fr.box, state_global->box);
818 copy_mat(state_global->box, state->box);
820 if (vsite && (Flags & MD_RERUN_VSITE))
822 if (DOMAINDECOMP(cr))
824 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
828 /* Following is necessary because the graph may get out of sync
829 * with the coordinates if we only have every N'th coordinate set
831 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
832 shift_self(graph, state->box, state->x);
834 construct_vsites(vsite, state->x, ir->delta_t, state->v,
835 top->idef.iparams, top->idef.il,
836 fr->ePBC, fr->bMolPBC, cr, state->box);
839 unshift_self(graph, state->box, state->x);
844 /* Stop Center of Mass motion */
845 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
849 /* for rerun MD always do Neighbour Searching */
850 bNS = (bFirstStep || ir->nstlist != 0);
855 /* Determine whether or not to do Neighbour Searching and LR */
856 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
858 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
859 (ir->nstlist == -1 && nlh.nabnsb > 0));
861 if (bNS && ir->nstlist == -1)
863 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
867 /* check whether we should stop because another simulation has
871 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
872 (multisim_nsteps != ir->nsteps) )
879 "Stopping simulation %d because another one has finished\n",
883 gs.sig[eglsCHKPT] = 1;
888 /* < 0 means stop at next step, > 0 means stop at next NS step */
889 if ( (gs.set[eglsSTOPCOND] < 0) ||
890 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
895 /* Determine whether or not to update the Born radii if doing GB */
896 bBornRadii = bFirstStep;
897 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
902 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
903 do_verbose = bVerbose &&
904 (step % stepout == 0 || bFirstStep || bLastStep);
906 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
914 bMasterState = FALSE;
915 /* Correct the new box if it is too skewed */
916 if (DYNAMIC_BOX(*ir))
918 if (correct_box(fplog, step, state->box, graph))
923 if (DOMAINDECOMP(cr) && bMasterState)
925 dd_collect_state(cr->dd, state, state_global);
929 if (DOMAINDECOMP(cr))
931 /* Repartition the domain decomposition */
932 wallcycle_start(wcycle, ewcDOMDEC);
933 dd_partition_system(fplog, step, cr,
934 bMasterState, nstglobalcomm,
935 state_global, top_global, ir,
936 state, &f, mdatoms, top, fr,
937 vsite, shellfc, constr,
939 do_verbose && !bPMETuneRunning);
940 wallcycle_stop(wcycle, ewcDOMDEC);
941 /* If using an iterative integrator, reallocate space to match the decomposition */
945 if (MASTER(cr) && do_log)
947 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
950 if (ir->efep != efepNO)
952 update_mdatoms(mdatoms, state->lambda[efptMASS]);
955 if ((bRerunMD && rerun_fr.bV) || bExchanged)
958 /* We need the kinetic energy at minus the half step for determining
959 * the full step kinetic energy and possibly for T-coupling.*/
960 /* This may not be quite working correctly yet . . . . */
961 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
962 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
963 constr, NULL, FALSE, state->box,
964 top_global, &bSumEkinhOld,
965 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
967 clear_mat(force_vir);
969 /* We write a checkpoint at this MD step when:
970 * either at an NS step when we signalled through gs,
971 * or at the last step (but not when we do not want confout),
972 * but never at the first step or with rerun.
974 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
975 (bLastStep && (Flags & MD_CONFOUT))) &&
976 step > ir->init_step && !bRerunMD);
979 gs.set[eglsCHKPT] = 0;
982 /* Determine the energy and pressure:
983 * at nstcalcenergy steps and at energy output steps (set below).
985 if (EI_VV(ir->eI) && (!bInitStep))
987 /* for vv, the first half of the integration actually corresponds
988 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
989 but the virial needs to be calculated on both the current step and the 'next' step. Future
990 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
992 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
993 bCalcVir = bCalcEner ||
994 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
998 bCalcEner = do_per_step(step, ir->nstcalcenergy);
999 bCalcVir = bCalcEner ||
1000 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1003 /* Do we need global communication ? */
1004 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1005 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1006 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1008 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1010 if (do_ene || do_log || bDoReplEx)
1017 /* these CGLO_ options remain the same throughout the iteration */
1018 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1019 (bGStat ? CGLO_GSTAT : 0)
1022 force_flags = (GMX_FORCE_STATECHANGED |
1023 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1024 GMX_FORCE_ALLFORCES |
1026 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1027 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1028 (bDoFEP ? GMX_FORCE_DHDL : 0)
1033 if (do_per_step(step, ir->nstcalclr))
1035 force_flags |= GMX_FORCE_DO_LR;
1041 /* Now is the time to relax the shells */
1042 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1043 ir, bNS, force_flags,
1046 state, f, force_vir, mdatoms,
1047 nrnb, wcycle, graph, groups,
1048 shellfc, fr, bBornRadii, t, mu_tot,
1050 mdoutf_get_fp_field(outf));
1060 /* The coordinates (x) are shifted (to get whole molecules)
1062 * This is parallellized as well, and does communication too.
1063 * Check comments in sim_util.c
1065 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1066 state->box, state->x, &state->hist,
1067 f, force_vir, mdatoms, enerd, fcd,
1068 state->lambda, graph,
1069 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1070 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1073 if (bVV && !bStartingFromCpt && !bRerunMD)
1074 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1076 if (ir->eI == eiVV && bInitStep)
1078 /* if using velocity verlet with full time step Ekin,
1079 * take the first half step only to compute the
1080 * virial for the first step. From there,
1081 * revert back to the initial coordinates
1082 * so that the input is actually the initial step.
1084 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1088 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1089 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1092 /* If we are using twin-range interactions where the long-range component
1093 * is only evaluated every nstcalclr>1 steps, we should do a special update
1094 * step to combine the long-range forces on these steps.
1095 * For nstcalclr=1 this is not done, since the forces would have been added
1096 * directly to the short-range forces already.
1098 * TODO Remove various aspects of VV+twin-range in master
1099 * branch, because VV integrators did not ever support
1100 * twin-range multiple time stepping with constraints.
1102 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1104 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1105 f, bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1106 ekind, M, upd, bInitStep, etrtVELOCITY1,
1107 cr, nrnb, constr, &top->idef);
1109 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1111 gmx_iterate_init(&iterate, TRUE);
1113 /* for iterations, we save these vectors, as we will be self-consistently iterating
1116 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1118 /* save the state */
1119 if (iterate.bIterationActive)
1121 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1124 bFirstIterate = TRUE;
1125 while (bFirstIterate || iterate.bIterationActive)
1127 if (iterate.bIterationActive)
1129 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1130 if (bFirstIterate && bTrotter)
1132 /* The first time through, we need a decent first estimate
1133 of veta(t+dt) to compute the constraints. Do
1134 this by computing the box volume part of the
1135 trotter integration at this time. Nothing else
1136 should be changed by this routine here. If
1137 !(first time), we start with the previous value
1140 veta_save = state->veta;
1141 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1142 vetanew = state->veta;
1143 state->veta = veta_save;
1148 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1150 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1151 state, fr->bMolPBC, graph, f,
1152 &top->idef, shake_vir,
1153 cr, nrnb, wcycle, upd, constr,
1154 TRUE, bCalcVir, vetanew);
1156 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1158 /* Correct the virial for multiple time stepping */
1159 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1164 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1170 /* Need to unshift here if a do_force has been
1171 called in the previous step */
1172 unshift_self(graph, state->box, state->x);
1175 /* if VV, compute the pressure and constraints */
1176 /* For VV2, we strictly only need this if using pressure
1177 * control, but we really would like to have accurate pressures
1179 * Think about ways around this in the future?
1180 * For now, keep this choice in comments.
1182 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1183 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1185 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1186 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1188 bSumEkinhOld = TRUE;
1190 /* for vv, the first half of the integration actually corresponds to the previous step.
1191 So we need information from the last step in the first half of the integration */
1192 if (bGStat || do_per_step(step-1, nstglobalcomm))
1194 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1195 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1196 constr, NULL, FALSE, state->box,
1197 top_global, &bSumEkinhOld,
1200 | (bTemp ? CGLO_TEMPERATURE : 0)
1201 | (bPres ? CGLO_PRESSURE : 0)
1202 | (bPres ? CGLO_CONSTRAINT : 0)
1203 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1204 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1207 /* explanation of above:
1208 a) We compute Ekin at the full time step
1209 if 1) we are using the AveVel Ekin, and it's not the
1210 initial step, or 2) if we are using AveEkin, but need the full
1211 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1212 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1213 EkinAveVel because it's needed for the pressure */
1215 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1220 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1221 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1228 /* We need the kinetic energy at minus the half step for determining
1229 * the full step kinetic energy and possibly for T-coupling.*/
1230 /* This may not be quite working correctly yet . . . . */
1231 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1232 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1233 constr, NULL, FALSE, state->box,
1234 top_global, &bSumEkinhOld,
1235 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1240 if (iterate.bIterationActive &&
1241 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1242 state->veta, &vetanew))
1246 bFirstIterate = FALSE;
1249 if (bTrotter && !bInitStep)
1251 copy_mat(shake_vir, state->svir_prev);
1252 copy_mat(force_vir, state->fvir_prev);
1253 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1255 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1256 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1257 enerd->term[F_EKIN] = trace(ekind->ekin);
1260 /* if it's the initial step, we performed this first step just to get the constraint virial */
1261 if (bInitStep && ir->eI == eiVV)
1263 copy_rvecn(cbuf, state->v, 0, state->natoms);
1267 /* MRS -- now done iterating -- compute the conserved quantity */
1270 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1273 last_ekin = enerd->term[F_EKIN];
1275 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1277 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1279 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1282 sum_dhdl(enerd, state->lambda, ir->fepvals);
1286 /* ######## END FIRST UPDATE STEP ############## */
1287 /* ######## If doing VV, we now have v(dt) ###### */
1290 /* perform extended ensemble sampling in lambda - we don't
1291 actually move to the new state before outputting
1292 statistics, but if performing simulated tempering, we
1293 do update the velocities and the tau_t. */
1295 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1296 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1297 copy_df_history(&state_global->dfhist, &state->dfhist);
1300 /* Now we have the energies and forces corresponding to the
1301 * coordinates at time t. We must output all of this before
1304 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1305 ir, state, state_global, top_global, fr,
1306 outf, mdebin, ekind, f, f_global,
1308 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1310 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1311 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1313 /* kludge -- virial is lost with restart for NPT control. Must restart */
1314 if (bStartingFromCpt && bVV)
1316 copy_mat(state->svir_prev, shake_vir);
1317 copy_mat(state->fvir_prev, force_vir);
1320 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1322 /* Check whether everything is still allright */
1323 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1324 #ifdef GMX_THREAD_MPI
1329 /* this is just make gs.sig compatible with the hack
1330 of sending signals around by MPI_Reduce with together with
1332 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1334 gs.sig[eglsSTOPCOND] = 1;
1336 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1338 gs.sig[eglsSTOPCOND] = -1;
1340 /* < 0 means stop at next step, > 0 means stop at next NS step */
1344 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1345 gmx_get_signal_name(),
1346 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1350 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1351 gmx_get_signal_name(),
1352 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1354 handled_stop_condition = (int)gmx_get_stop_condition();
1356 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1357 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1358 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1360 /* Signal to terminate the run */
1361 gs.sig[eglsSTOPCOND] = 1;
1364 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1366 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1369 if (bResetCountersHalfMaxH && MASTER(cr) &&
1370 elapsed_time > max_hours*60.0*60.0*0.495)
1372 gs.sig[eglsRESETCOUNTERS] = 1;
1375 if (ir->nstlist == -1 && !bRerunMD)
1377 /* When bGStatEveryStep=FALSE, global_stat is only called
1378 * when we check the atom displacements, not at NS steps.
1379 * This means that also the bonded interaction count check is not
1380 * performed immediately after NS. Therefore a few MD steps could
1381 * be performed with missing interactions.
1382 * But wrong energies are never written to file,
1383 * since energies are only written after global_stat
1386 if (step >= nlh.step_nscheck)
1388 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1389 nlh.scale_tot, state->x);
1393 /* This is not necessarily true,
1394 * but step_nscheck is determined quite conservatively.
1400 /* In parallel we only have to check for checkpointing in steps
1401 * where we do global communication,
1402 * otherwise the other nodes don't know.
1404 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1407 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1408 gs.set[eglsCHKPT] == 0)
1410 gs.sig[eglsCHKPT] = 1;
1413 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1418 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1420 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1422 gmx_bool bIfRandomize;
1423 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1424 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1425 if (constr && bIfRandomize)
1427 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1428 state, fr->bMolPBC, graph, f,
1429 &top->idef, tmp_vir,
1430 cr, nrnb, wcycle, upd, constr,
1431 TRUE, bCalcVir, vetanew);
1436 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1438 gmx_iterate_init(&iterate, TRUE);
1439 /* for iterations, we save these vectors, as we will be redoing the calculations */
1440 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1443 bFirstIterate = TRUE;
1444 while (bFirstIterate || iterate.bIterationActive)
1446 /* We now restore these vectors to redo the calculation with improved extended variables */
1447 if (iterate.bIterationActive)
1449 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1452 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1453 so scroll down for that logic */
1455 /* ######### START SECOND UPDATE STEP ################# */
1456 /* Box is changed in update() when we do pressure coupling,
1457 * but we should still use the old box for energy corrections and when
1458 * writing it to the energy file, so it matches the trajectory files for
1459 * the same timestep above. Make a copy in a separate array.
1461 copy_mat(state->box, lastbox);
1466 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1468 wallcycle_start(wcycle, ewcUPDATE);
1469 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1472 if (iterate.bIterationActive)
1480 /* we use a new value of scalevir to converge the iterations faster */
1481 scalevir = tracevir/trace(shake_vir);
1483 msmul(shake_vir, scalevir, shake_vir);
1484 m_add(force_vir, shake_vir, total_vir);
1485 clear_mat(shake_vir);
1487 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1488 /* We can only do Berendsen coupling after we have summed
1489 * the kinetic energy or virial. Since the happens
1490 * in global_state after update, we should only do it at
1491 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1496 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1497 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1502 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1504 /* velocity half-step update */
1505 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1506 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1507 ekind, M, upd, FALSE, etrtVELOCITY2,
1508 cr, nrnb, constr, &top->idef);
1511 /* Above, initialize just copies ekinh into ekin,
1512 * it doesn't copy position (for VV),
1513 * and entire integrator for MD.
1516 if (ir->eI == eiVVAK)
1518 copy_rvecn(state->x, cbuf, 0, state->natoms);
1520 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1522 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1523 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1524 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1525 wallcycle_stop(wcycle, ewcUPDATE);
1527 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1528 fr->bMolPBC, graph, f,
1529 &top->idef, shake_vir,
1530 cr, nrnb, wcycle, upd, constr,
1531 FALSE, bCalcVir, state->veta);
1533 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1535 /* Correct the virial for multiple time stepping */
1536 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1539 if (ir->eI == eiVVAK)
1541 /* erase F_EKIN and F_TEMP here? */
1542 /* just compute the kinetic energy at the half step to perform a trotter step */
1543 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1544 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1545 constr, NULL, FALSE, lastbox,
1546 top_global, &bSumEkinhOld,
1547 cglo_flags | CGLO_TEMPERATURE
1549 wallcycle_start(wcycle, ewcUPDATE);
1550 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1551 /* now we know the scaling, we can compute the positions again again */
1552 copy_rvecn(cbuf, state->x, 0, state->natoms);
1554 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1556 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1557 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1558 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1559 wallcycle_stop(wcycle, ewcUPDATE);
1561 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1562 /* are the small terms in the shake_vir here due
1563 * to numerical errors, or are they important
1564 * physically? I'm thinking they are just errors, but not completely sure.
1565 * For now, will call without actually constraining, constr=NULL*/
1566 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1567 state, fr->bMolPBC, graph, f,
1568 &top->idef, tmp_vir,
1569 cr, nrnb, wcycle, upd, NULL,
1575 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1578 if (fr->bSepDVDL && fplog && do_log)
1580 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1584 /* this factor or 2 correction is necessary
1585 because half of the constraint force is removed
1586 in the vv step, so we have to double it. See
1587 the Redmine issue #1255. It is not yet clear
1588 if the factor of 2 is exact, or just a very
1589 good approximation, and this will be
1590 investigated. The next step is to see if this
1591 can be done adding a dhdl contribution from the
1592 rattle step, but this is somewhat more
1593 complicated with the current code. Will be
1594 investigated, hopefully for 4.6.3. However,
1595 this current solution is much better than
1596 having it completely wrong.
1598 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1602 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1607 /* Need to unshift here */
1608 unshift_self(graph, state->box, state->x);
1613 wallcycle_start(wcycle, ewcVSITECONSTR);
1616 shift_self(graph, state->box, state->x);
1618 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1619 top->idef.iparams, top->idef.il,
1620 fr->ePBC, fr->bMolPBC, cr, state->box);
1624 unshift_self(graph, state->box, state->x);
1626 wallcycle_stop(wcycle, ewcVSITECONSTR);
1629 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1630 /* With Leap-Frog we can skip compute_globals at
1631 * non-communication steps, but we need to calculate
1632 * the kinetic energy one step before communication.
1634 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1636 if (ir->nstlist == -1 && bFirstIterate)
1638 gs.sig[eglsNABNSB] = nlh.nabnsb;
1640 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1641 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1643 bFirstIterate ? &gs : NULL,
1644 (step_rel % gs.nstms == 0) &&
1645 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1647 top_global, &bSumEkinhOld,
1649 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1650 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1651 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1652 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1653 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1654 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1657 if (ir->nstlist == -1 && bFirstIterate)
1659 nlh.nabnsb = gs.set[eglsNABNSB];
1660 gs.set[eglsNABNSB] = 0;
1663 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1664 /* ############# END CALC EKIN AND PRESSURE ################# */
1666 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1667 the virial that should probably be addressed eventually. state->veta has better properies,
1668 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1669 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1671 if (iterate.bIterationActive &&
1672 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1673 trace(shake_vir), &tracevir))
1677 bFirstIterate = FALSE;
1680 if (!bVV || bRerunMD)
1682 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1683 sum_dhdl(enerd, state->lambda, ir->fepvals);
1685 update_box(fplog, step, ir, mdatoms, state, f,
1686 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1688 /* ################# END UPDATE STEP 2 ################# */
1689 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1691 /* The coordinates (x) were unshifted in update */
1694 /* We will not sum ekinh_old,
1695 * so signal that we still have to do it.
1697 bSumEkinhOld = TRUE;
1700 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1702 /* use the directly determined last velocity, not actually the averaged half steps */
1703 if (bTrotter && ir->eI == eiVV)
1705 enerd->term[F_EKIN] = last_ekin;
1707 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1711 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1715 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1717 /* ######### END PREPARING EDR OUTPUT ########### */
1722 gmx_bool do_dr, do_or;
1724 if (fplog && do_log && bDoExpanded)
1726 /* only needed if doing expanded ensemble */
1727 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1728 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1730 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1734 upd_mdebin(mdebin, bDoDHDL, TRUE,
1735 t, mdatoms->tmass, enerd, state,
1736 ir->fepvals, ir->expandedvals, lastbox,
1737 shake_vir, force_vir, total_vir, pres,
1738 ekind, mu_tot, constr);
1742 upd_mdebin_step(mdebin);
1745 do_dr = do_per_step(step, ir->nstdisreout);
1746 do_or = do_per_step(step, ir->nstorireout);
1748 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1750 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1752 if (ir->ePull != epullNO)
1754 pull_print_output(ir->pull, step, t);
1757 if (do_per_step(step, ir->nstlog))
1759 if (fflush(fplog) != 0)
1761 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1767 /* Have to do this part _after_ outputting the logfile and the edr file */
1768 /* Gets written into the state at the beginning of next loop*/
1769 state->fep_state = lamnew;
1771 /* Print the remaining wall clock time for the run */
1772 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1776 fprintf(stderr, "\n");
1778 print_time(stderr, walltime_accounting, step, ir, cr);
1781 /* Ion/water position swapping.
1782 * Not done in last step since trajectory writing happens before this call
1783 * in the MD loop and exchanges would be lost anyway. */
1784 bNeedRepartition = FALSE;
1785 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1786 do_per_step(step, ir->swap->nstswap))
1788 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1789 bRerunMD ? rerun_fr.x : state->x,
1790 bRerunMD ? rerun_fr.box : state->box,
1791 top_global, MASTER(cr) && bVerbose, bRerunMD);
1793 if (bNeedRepartition && DOMAINDECOMP(cr))
1795 dd_collect_state(cr->dd, state, state_global);
1799 /* Replica exchange */
1803 bExchanged = replica_exchange(fplog, cr, repl_ex,
1804 state_global, enerd,
1808 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1810 dd_partition_system(fplog, step, cr, TRUE, 1,
1811 state_global, top_global, ir,
1812 state, &f, mdatoms, top, fr,
1813 vsite, shellfc, constr,
1814 nrnb, wcycle, FALSE);
1819 bStartingFromCpt = FALSE;
1821 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1822 /* With all integrators, except VV, we need to retain the pressure
1823 * at the current step for coupling at the next step.
1825 if ((state->flags & (1<<estPRES_PREV)) &&
1827 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1829 /* Store the pressure in t_state for pressure coupling
1830 * at the next MD step.
1832 copy_mat(pres, state->pres_prev);
1835 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1837 if ( (membed != NULL) && (!bLastStep) )
1839 rescale_membed(step_rel, membed, state_global->x);
1846 /* read next frame from input trajectory */
1847 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1852 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1856 if (!bRerunMD || !rerun_fr.bStep)
1858 /* increase the MD step number */
1863 cycles = wallcycle_stop(wcycle, ewcSTEP);
1864 if (DOMAINDECOMP(cr) && wcycle)
1866 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1869 if (bPMETuneRunning || bPMETuneTry)
1871 /* PME grid + cut-off optimization with GPUs or PME nodes */
1873 /* Count the total cycles over the last steps */
1874 cycles_pmes += cycles;
1876 /* We can only switch cut-off at NS steps */
1877 if (step % ir->nstlist == 0)
1879 /* PME grid + cut-off optimization with GPUs or PME nodes */
1882 if (DDMASTER(cr->dd))
1884 /* PME node load is too high, start tuning */
1885 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1887 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1889 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1891 bPMETuneTry = FALSE;
1894 if (bPMETuneRunning)
1896 /* init_step might not be a multiple of nstlist,
1897 * but the first cycle is always skipped anyhow.
1900 pme_load_balance(pme_loadbal, cr,
1901 (bVerbose && MASTER(cr)) ? stderr : NULL,
1903 ir, state, cycles_pmes,
1904 fr->ic, fr->nbv, &fr->pmedata,
1907 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1908 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1909 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1910 fr->rlist = fr->ic->rlist;
1911 fr->rlistlong = fr->ic->rlistlong;
1912 fr->rcoulomb = fr->ic->rcoulomb;
1913 fr->rvdw = fr->ic->rvdw;
1915 if (ir->eDispCorr != edispcNO)
1917 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1924 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1925 gs.set[eglsRESETCOUNTERS] != 0)
1927 /* Reset all the counters related to performance over the run */
1928 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1929 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1930 wcycle_set_reset_counters(wcycle, -1);
1931 if (!(cr->duty & DUTY_PME))
1933 /* Tell our PME node to reset its counters */
1934 gmx_pme_send_resetcounters(cr, step);
1936 /* Correct max_hours for the elapsed time */
1937 max_hours -= elapsed_time/(60.0*60.0);
1938 bResetCountersHalfMaxH = FALSE;
1939 gs.set[eglsRESETCOUNTERS] = 0;
1942 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1943 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1946 /* End of main MD loop */
1949 /* Stop measuring walltime */
1950 walltime_accounting_end(walltime_accounting);
1952 if (bRerunMD && MASTER(cr))
1957 if (!(cr->duty & DUTY_PME))
1959 /* Tell the PME only node to finish */
1960 gmx_pme_send_finish(cr);
1965 if (ir->nstcalcenergy > 0 && !bRerunMD)
1967 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1968 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1975 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1977 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)));
1978 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1981 if (pme_loadbal != NULL)
1983 pme_loadbal_done(pme_loadbal, cr, fplog,
1984 fr->nbv != NULL && fr->nbv->bUseGPU);
1987 if (shellfc && fplog)
1989 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1990 (nconverged*100.0)/step_rel);
1991 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1995 if (repl_ex_nst > 0 && MASTER(cr))
1997 print_replica_exchange_statistics(fplog, repl_ex);
2000 /* IMD cleanup, if bIMD is TRUE. */
2001 IMD_finalize(ir->bIMD, ir->imd);
2003 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);