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"
54 #include "md_support.h"
55 #include "md_logging.h"
67 #include "domdec_network.h"
68 #include "gromacs/gmxlib/topsort.h"
72 #include "gromacs/gmxpreprocess/compute_io.h"
73 #include "checkpoint.h"
74 #include "mtop_util.h"
75 #include "sighandler.h"
77 #include "gromacs/utility/cstringutil.h"
78 #include "pme_loadbal.h"
81 #include "types/nlistheuristics.h"
82 #include "types/iteratedconstraints.h"
83 #include "nbnxn_cuda_data_mgmt.h"
85 #include "gromacs/fileio/confio.h"
86 #include "gromacs/fileio/trajectory_writing.h"
87 #include "gromacs/fileio/trnio.h"
88 #include "gromacs/fileio/trxio.h"
89 #include "gromacs/fileio/xtcio.h"
90 #include "gromacs/imd/imd.h"
91 #include "gromacs/pbcutil/mshift.h"
92 #include "gromacs/pbcutil/pbc.h"
93 #include "gromacs/pulling/pull.h"
94 #include "gromacs/swap/swapcoords.h"
95 #include "gromacs/timing/wallcycle.h"
96 #include "gromacs/timing/walltime_accounting.h"
97 #include "gromacs/utility/gmxmpi.h"
98 #include "gromacs/utility/smalloc.h"
101 #include "corewrap.h"
104 static void reset_all_counters(FILE *fplog, t_commrec *cr,
106 gmx_int64_t *step_rel, t_inputrec *ir,
107 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
108 gmx_walltime_accounting_t walltime_accounting,
109 nbnxn_cuda_ptr_t cu_nbv)
111 char sbuf[STEPSTRSIZE];
113 /* Reset all the counters related to performance over the run */
114 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
115 gmx_step_str(step, sbuf));
119 nbnxn_cuda_reset_timings(cu_nbv);
122 wallcycle_stop(wcycle, ewcRUN);
123 wallcycle_reset_all(wcycle);
124 if (DOMAINDECOMP(cr))
126 reset_dd_statistics_counters(cr->dd);
129 ir->init_step += *step_rel;
130 ir->nsteps -= *step_rel;
132 wallcycle_start(wcycle, ewcRUN);
133 walltime_accounting_start(walltime_accounting);
134 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
137 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
138 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
140 gmx_vsite_t *vsite, gmx_constr_t constr,
141 int stepout, t_inputrec *ir,
142 gmx_mtop_t *top_global,
144 t_state *state_global,
146 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
147 gmx_edsam_t ed, t_forcerec *fr,
148 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
149 real cpt_period, real max_hours,
150 const char gmx_unused *deviceOptions,
153 gmx_walltime_accounting_t walltime_accounting)
155 gmx_mdoutf_t outf = NULL;
156 gmx_int64_t step, step_rel;
158 double t, t0, lam0[efptNR];
159 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
160 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
161 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
162 bBornRadii, bStartingFromCpt;
163 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
164 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
165 bForceUpdate = FALSE, bCPT;
166 gmx_bool bMasterState;
167 int force_flags, cglo_flags;
168 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
173 t_state *bufstate = NULL;
174 matrix *scale_tot, pcoupl_mu, M, ebox;
177 gmx_repl_ex_t repl_ex = NULL;
180 t_mdebin *mdebin = NULL;
181 t_state *state = NULL;
182 rvec *f_global = NULL;
183 gmx_enerdata_t *enerd;
185 gmx_global_stat_t gstat;
186 gmx_update_t upd = NULL;
187 t_graph *graph = NULL;
189 gmx_groups_t *groups;
190 gmx_ekindata_t *ekind, *ekind_save;
191 gmx_shellfc_t shellfc;
192 int count, nconverged = 0;
195 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged, bNeedRepartition;
197 gmx_bool bResetCountersHalfMaxH = FALSE;
198 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
199 gmx_bool bUpdateDoLR;
203 real veta_save, scalevir, tracevir;
209 real saved_conserved_quantity = 0;
214 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
215 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
216 gmx_iterate_t iterate;
217 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
218 simulation stops. If equal to zero, don't
219 communicate any more between multisims.*/
220 /* PME load balancing data for GPU kernels */
221 pme_load_balancing_t pme_loadbal = NULL;
223 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
226 gmx_bool bIMDstep = FALSE;
229 /* Temporary addition for FAHCORE checkpointing */
233 /* Check for special mdrun options */
234 bRerunMD = (Flags & MD_RERUN);
235 bAppend = (Flags & MD_APPENDFILES);
236 if (Flags & MD_RESETCOUNTERSHALFWAY)
240 /* Signal to reset the counters half the simulation steps. */
241 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
243 /* Signal to reset the counters halfway the simulation time. */
244 bResetCountersHalfMaxH = (max_hours > 0);
247 /* md-vv uses averaged full step velocities for T-control
248 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
249 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
251 if (bVV) /* to store the initial velocities while computing virial */
253 snew(cbuf, top_global->natoms);
255 /* all the iteratative cases - only if there are constraints */
256 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
257 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
258 false in this step. The correct value, true or false,
259 is set at each step, as it depends on the frequency of temperature
260 and pressure control.*/
261 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
265 /* Since we don't know if the frames read are related in any way,
266 * rebuild the neighborlist at every step.
269 ir->nstcalcenergy = 1;
273 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
275 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
276 bGStatEveryStep = (nstglobalcomm == 1);
278 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
281 "To reduce the energy communication with nstlist = -1\n"
282 "the neighbor list validity should not be checked at every step,\n"
283 "this means that exact integration is not guaranteed.\n"
284 "The neighbor list validity is checked after:\n"
285 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
286 "In most cases this will result in exact integration.\n"
287 "This reduces the energy communication by a factor of 2 to 3.\n"
288 "If you want less energy communication, set nstlist > 3.\n\n");
293 ir->nstxout_compressed = 0;
295 groups = &top_global->groups;
298 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
299 &(state_global->fep_state), lam0,
300 nrnb, top_global, &upd,
301 nfile, fnm, &outf, &mdebin,
302 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
304 clear_mat(total_vir);
306 /* Energy terms and groups */
308 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
310 if (DOMAINDECOMP(cr))
316 snew(f, top_global->natoms);
319 /* Kinetic energy data */
321 init_ekindata(fplog, top_global, &(ir->opts), ekind);
322 /* needed for iteration of constraints */
324 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
325 /* Copy the cos acceleration to the groups struct */
326 ekind->cosacc.cos_accel = ir->cos_accel;
328 gstat = global_stat_init(ir);
331 /* Check for polarizable models and flexible constraints */
332 shellfc = init_shell_flexcon(fplog, fr->cutoff_scheme == ecutsVERLET,
333 top_global, n_flexible_constraints(constr),
334 (ir->bContinuation ||
335 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
336 NULL : state_global->x);
338 if (shellfc && ir->eI == eiNM)
340 /* Currently shells don't work with Normal Modes */
341 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");
344 if (vsite && ir->eI == eiNM)
346 /* Currently virtual sites don't work with Normal Modes */
347 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");
352 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
353 set_deform_reference_box(upd,
354 deform_init_init_step_tpx,
355 deform_init_box_tpx);
356 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
360 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
361 if ((io > 2000) && MASTER(cr))
364 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
369 if (DOMAINDECOMP(cr))
371 top = dd_init_local_top(top_global);
374 dd_init_local_state(cr->dd, state_global, state);
376 if (DDMASTER(cr->dd) && ir->nstfout)
378 snew(f_global, state_global->natoms);
383 top = gmx_mtop_generate_local_top(top_global, ir);
385 forcerec_set_excl_load(fr, top);
387 state = serial_init_local_state(state_global);
390 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
394 set_vsite_top(vsite, top, mdatoms, cr);
397 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
399 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
404 make_local_shells(cr, mdatoms, shellfc);
407 setup_bonded_threading(fr, &top->idef);
410 /* Set up interactive MD (IMD) */
411 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
412 nfile, fnm, oenv, imdport, Flags);
414 if (DOMAINDECOMP(cr))
416 /* Distribute the charge groups over the nodes from the master node */
417 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
418 state_global, top_global, ir,
419 state, &f, mdatoms, top, fr,
420 vsite, shellfc, constr,
421 nrnb, wcycle, FALSE);
425 update_mdatoms(mdatoms, state->lambda[efptMASS]);
427 if (opt2bSet("-cpi", nfile, fnm))
429 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
433 bStateFromCP = FALSE;
438 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
445 /* Update mdebin with energy history if appending to output files */
446 if (Flags & MD_APPENDFILES)
448 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
452 /* We might have read an energy history from checkpoint,
453 * free the allocated memory and reset the counts.
455 done_energyhistory(&state_global->enerhist);
456 init_energyhistory(&state_global->enerhist);
459 /* Set the initial energy history in state by updating once */
460 update_energyhistory(&state_global->enerhist, mdebin);
463 /* Initialize constraints */
464 if (constr && !DOMAINDECOMP(cr))
466 set_constraints(constr, top, ir, mdatoms, cr);
471 /* We need to be sure replica exchange can only occur
472 * when the energies are current */
473 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
474 "repl_ex_nst", &repl_ex_nst);
475 /* This check needs to happen before inter-simulation
476 * signals are initialized, too */
478 if (repl_ex_nst > 0 && MASTER(cr))
480 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
481 repl_ex_nst, repl_ex_nex, repl_ex_seed);
484 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
485 * PME tuning is not supported with PME only for LJ and not for Coulomb.
487 if ((Flags & MD_TUNEPME) &&
488 EEL_PME(fr->eeltype) &&
489 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
492 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
494 if (cr->duty & DUTY_PME)
496 /* Start tuning right away, as we can't measure the load */
497 bPMETuneRunning = TRUE;
501 /* Separate PME nodes, we can measure the PP/PME load balance */
506 if (!ir->bContinuation && !bRerunMD)
508 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
510 /* Set the velocities of frozen particles to zero */
511 for (i = 0; i < mdatoms->homenr; i++)
513 for (m = 0; m < DIM; m++)
515 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
525 /* Constrain the initial coordinates and velocities */
526 do_constrain_first(fplog, constr, ir, mdatoms, state,
531 /* Construct the virtual sites for the initial configuration */
532 construct_vsites(vsite, state->x, ir->delta_t, NULL,
533 top->idef.iparams, top->idef.il,
534 fr->ePBC, fr->bMolPBC, cr, state->box);
540 /* set free energy calculation frequency as the minimum
541 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
542 nstfep = ir->fepvals->nstdhdl;
545 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
549 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
552 /* I'm assuming we need global communication the first time! MRS */
553 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
554 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
555 | (bVV ? CGLO_PRESSURE : 0)
556 | (bVV ? CGLO_CONSTRAINT : 0)
557 | (bRerunMD ? CGLO_RERUNMD : 0)
558 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
560 bSumEkinhOld = FALSE;
561 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
562 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
563 constr, NULL, FALSE, state->box,
564 top_global, &bSumEkinhOld, cglo_flags);
565 if (ir->eI == eiVVAK)
567 /* a second call to get the half step temperature initialized as well */
568 /* we do the same call as above, but turn the pressure off -- internally to
569 compute_globals, this is recognized as a velocity verlet half-step
570 kinetic energy calculation. This minimized excess variables, but
571 perhaps loses some logic?*/
573 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
574 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
575 constr, NULL, FALSE, state->box,
576 top_global, &bSumEkinhOld,
577 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
580 /* Calculate the initial half step temperature, and save the ekinh_old */
581 if (!(Flags & MD_STARTFROMCPT))
583 for (i = 0; (i < ir->opts.ngtc); i++)
585 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
590 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
591 and there is no previous step */
594 /* if using an iterative algorithm, we need to create a working directory for the state. */
597 bufstate = init_bufstate(state);
600 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
601 temperature control */
602 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
606 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
609 "RMS relative constraint deviation after constraining: %.2e\n",
610 constr_rmsd(constr, FALSE));
612 if (EI_STATE_VELOCITY(ir->eI))
614 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
618 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
619 " input trajectory '%s'\n\n",
620 *(top_global->name), opt2fn("-rerun", nfile, fnm));
623 fprintf(stderr, "Calculated time to finish depends on nsteps from "
624 "run input file,\nwhich may not correspond to the time "
625 "needed to process input trajectory.\n\n");
631 fprintf(stderr, "starting mdrun '%s'\n",
632 *(top_global->name));
635 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
639 sprintf(tbuf, "%s", "infinite");
641 if (ir->init_step > 0)
643 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
644 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
645 gmx_step_str(ir->init_step, sbuf2),
646 ir->init_step*ir->delta_t);
650 fprintf(stderr, "%s steps, %s ps.\n",
651 gmx_step_str(ir->nsteps, sbuf), tbuf);
654 fprintf(fplog, "\n");
657 walltime_accounting_start(walltime_accounting);
658 wallcycle_start(wcycle, ewcRUN);
659 print_start(fplog, cr, walltime_accounting, "mdrun");
661 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
663 chkpt_ret = fcCheckPointParallel( cr->nodeid,
667 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
672 /***********************************************************
676 ************************************************************/
678 /* if rerunMD then read coordinates and velocities from input trajectory */
681 if (getenv("GMX_FORCE_UPDATE"))
689 bNotLastFrame = read_first_frame(oenv, &status,
690 opt2fn("-rerun", nfile, fnm),
691 &rerun_fr, TRX_NEED_X | TRX_READ_V);
692 if (rerun_fr.natoms != top_global->natoms)
695 "Number of atoms in trajectory (%d) does not match the "
696 "run input file (%d)\n",
697 rerun_fr.natoms, top_global->natoms);
699 if (ir->ePBC != epbcNONE)
703 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);
705 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
707 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
714 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
717 if (ir->ePBC != epbcNONE)
719 /* Set the shift vectors.
720 * Necessary here when have a static box different from the tpr box.
722 calc_shifts(rerun_fr.box, fr->shift_vec);
726 /* loop over MD steps or if rerunMD to end of input trajectory */
728 /* Skip the first Nose-Hoover integration when we get the state from tpx */
729 bStateFromTPX = !bStateFromCP;
730 bInitStep = bFirstStep && (bStateFromTPX || bVV);
731 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
733 bSumEkinhOld = FALSE;
735 bNeedRepartition = FALSE;
737 init_global_signals(&gs, cr, ir, repl_ex_nst);
739 step = ir->init_step;
742 if (ir->nstlist == -1)
744 init_nlistheuristics(&nlh, bGStatEveryStep, step);
747 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
749 /* check how many steps are left in other sims */
750 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
754 /* and stop now if we should */
755 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
756 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
757 while (!bLastStep || (bRerunMD && bNotLastFrame))
760 wallcycle_start(wcycle, ewcSTEP);
766 step = rerun_fr.step;
767 step_rel = step - ir->init_step;
780 bLastStep = (step_rel == ir->nsteps);
781 t = t0 + step*ir->delta_t;
784 if (ir->efep != efepNO || ir->bSimTemp)
786 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
787 requiring different logic. */
789 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
790 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
791 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
792 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
793 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
798 update_annealing_target_temp(&(ir->opts), t);
803 if (!DOMAINDECOMP(cr) || MASTER(cr))
805 for (i = 0; i < state_global->natoms; i++)
807 copy_rvec(rerun_fr.x[i], state_global->x[i]);
811 for (i = 0; i < state_global->natoms; i++)
813 copy_rvec(rerun_fr.v[i], state_global->v[i]);
818 for (i = 0; i < state_global->natoms; i++)
820 clear_rvec(state_global->v[i]);
824 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
825 " Ekin, temperature and pressure are incorrect,\n"
826 " the virial will be incorrect when constraints are present.\n"
828 bRerunWarnNoV = FALSE;
832 copy_mat(rerun_fr.box, state_global->box);
833 copy_mat(state_global->box, state->box);
835 if (vsite && (Flags & MD_RERUN_VSITE))
837 if (DOMAINDECOMP(cr))
839 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
843 /* Following is necessary because the graph may get out of sync
844 * with the coordinates if we only have every N'th coordinate set
846 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
847 shift_self(graph, state->box, state->x);
849 construct_vsites(vsite, state->x, ir->delta_t, state->v,
850 top->idef.iparams, top->idef.il,
851 fr->ePBC, fr->bMolPBC, cr, state->box);
854 unshift_self(graph, state->box, state->x);
859 /* Stop Center of Mass motion */
860 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
864 /* for rerun MD always do Neighbour Searching */
865 bNS = (bFirstStep || ir->nstlist != 0);
870 /* Determine whether or not to do Neighbour Searching and LR */
871 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
873 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
874 (ir->nstlist == -1 && nlh.nabnsb > 0));
876 if (bNS && ir->nstlist == -1)
878 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
882 /* check whether we should stop because another simulation has
886 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
887 (multisim_nsteps != ir->nsteps) )
894 "Stopping simulation %d because another one has finished\n",
898 gs.sig[eglsCHKPT] = 1;
903 /* < 0 means stop at next step, > 0 means stop at next NS step */
904 if ( (gs.set[eglsSTOPCOND] < 0) ||
905 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
910 /* Determine whether or not to update the Born radii if doing GB */
911 bBornRadii = bFirstStep;
912 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
917 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
918 do_verbose = bVerbose &&
919 (step % stepout == 0 || bFirstStep || bLastStep);
921 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
929 bMasterState = FALSE;
930 /* Correct the new box if it is too skewed */
931 if (DYNAMIC_BOX(*ir))
933 if (correct_box(fplog, step, state->box, graph))
938 if (DOMAINDECOMP(cr) && bMasterState)
940 dd_collect_state(cr->dd, state, state_global);
944 if (DOMAINDECOMP(cr))
946 /* Repartition the domain decomposition */
947 wallcycle_start(wcycle, ewcDOMDEC);
948 dd_partition_system(fplog, step, cr,
949 bMasterState, nstglobalcomm,
950 state_global, top_global, ir,
951 state, &f, mdatoms, top, fr,
952 vsite, shellfc, constr,
954 do_verbose && !bPMETuneRunning);
955 wallcycle_stop(wcycle, ewcDOMDEC);
956 /* If using an iterative integrator, reallocate space to match the decomposition */
960 if (MASTER(cr) && do_log)
962 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
965 if (ir->efep != efepNO)
967 update_mdatoms(mdatoms, state->lambda[efptMASS]);
970 if ((bRerunMD && rerun_fr.bV) || bExchanged)
973 /* We need the kinetic energy at minus the half step for determining
974 * the full step kinetic energy and possibly for T-coupling.*/
975 /* This may not be quite working correctly yet . . . . */
976 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
977 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
978 constr, NULL, FALSE, state->box,
979 top_global, &bSumEkinhOld,
980 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
982 clear_mat(force_vir);
984 /* We write a checkpoint at this MD step when:
985 * either at an NS step when we signalled through gs,
986 * or at the last step (but not when we do not want confout),
987 * but never at the first step or with rerun.
989 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
990 (bLastStep && (Flags & MD_CONFOUT))) &&
991 step > ir->init_step && !bRerunMD);
994 gs.set[eglsCHKPT] = 0;
997 /* Determine the energy and pressure:
998 * at nstcalcenergy steps and at energy output steps (set below).
1000 if (EI_VV(ir->eI) && (!bInitStep))
1002 /* for vv, the first half of the integration actually corresponds
1003 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1004 but the virial needs to be calculated on both the current step and the 'next' step. Future
1005 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1007 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1008 bCalcVir = bCalcEner ||
1009 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1013 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1014 bCalcVir = bCalcEner ||
1015 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1018 /* Do we need global communication ? */
1019 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1020 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1021 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1023 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1025 if (do_ene || do_log)
1032 /* these CGLO_ options remain the same throughout the iteration */
1033 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1034 (bGStat ? CGLO_GSTAT : 0)
1037 force_flags = (GMX_FORCE_STATECHANGED |
1038 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1039 GMX_FORCE_ALLFORCES |
1041 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1042 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1043 (bDoFEP ? GMX_FORCE_DHDL : 0)
1048 if (do_per_step(step, ir->nstcalclr))
1050 force_flags |= GMX_FORCE_DO_LR;
1056 /* Now is the time to relax the shells */
1057 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1058 ir, bNS, force_flags,
1061 state, f, force_vir, mdatoms,
1062 nrnb, wcycle, graph, groups,
1063 shellfc, fr, bBornRadii, t, mu_tot,
1065 mdoutf_get_fp_field(outf));
1075 /* The coordinates (x) are shifted (to get whole molecules)
1077 * This is parallellized as well, and does communication too.
1078 * Check comments in sim_util.c
1080 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1081 state->box, state->x, &state->hist,
1082 f, force_vir, mdatoms, enerd, fcd,
1083 state->lambda, graph,
1084 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1085 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1088 if (bVV && !bStartingFromCpt && !bRerunMD)
1089 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1091 if (ir->eI == eiVV && bInitStep)
1093 /* if using velocity verlet with full time step Ekin,
1094 * take the first half step only to compute the
1095 * virial for the first step. From there,
1096 * revert back to the initial coordinates
1097 * so that the input is actually the initial step.
1099 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1103 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1104 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1107 /* If we are using twin-range interactions where the long-range component
1108 * is only evaluated every nstcalclr>1 steps, we should do a special update
1109 * step to combine the long-range forces on these steps.
1110 * For nstcalclr=1 this is not done, since the forces would have been added
1111 * directly to the short-range forces already.
1113 * TODO Remove various aspects of VV+twin-range in master
1114 * branch, because VV integrators did not ever support
1115 * twin-range multiple time stepping with constraints.
1117 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1119 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1120 f, bUpdateDoLR, fr->f_twin, fcd,
1121 ekind, M, upd, bInitStep, etrtVELOCITY1,
1122 cr, nrnb, constr, &top->idef);
1124 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1126 gmx_iterate_init(&iterate, TRUE);
1128 /* for iterations, we save these vectors, as we will be self-consistently iterating
1131 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1133 /* save the state */
1134 if (iterate.bIterationActive)
1136 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1139 bFirstIterate = TRUE;
1140 while (bFirstIterate || iterate.bIterationActive)
1142 if (iterate.bIterationActive)
1144 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1145 if (bFirstIterate && bTrotter)
1147 /* The first time through, we need a decent first estimate
1148 of veta(t+dt) to compute the constraints. Do
1149 this by computing the box volume part of the
1150 trotter integration at this time. Nothing else
1151 should be changed by this routine here. If
1152 !(first time), we start with the previous value
1155 veta_save = state->veta;
1156 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1157 vetanew = state->veta;
1158 state->veta = veta_save;
1163 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1165 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1166 state, fr->bMolPBC, graph, f,
1167 &top->idef, shake_vir,
1168 cr, nrnb, wcycle, upd, constr,
1169 TRUE, bCalcVir, vetanew);
1173 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1179 /* Need to unshift here if a do_force has been
1180 called in the previous step */
1181 unshift_self(graph, state->box, state->x);
1184 /* if VV, compute the pressure and constraints */
1185 /* For VV2, we strictly only need this if using pressure
1186 * control, but we really would like to have accurate pressures
1188 * Think about ways around this in the future?
1189 * For now, keep this choice in comments.
1191 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1192 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1194 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1195 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1197 bSumEkinhOld = TRUE;
1199 /* for vv, the first half of the integration actually corresponds to the previous step.
1200 So we need information from the last step in the first half of the integration */
1201 if (bGStat || do_per_step(step-1, nstglobalcomm))
1203 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1204 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1205 constr, NULL, FALSE, state->box,
1206 top_global, &bSumEkinhOld,
1209 | (bTemp ? CGLO_TEMPERATURE : 0)
1210 | (bPres ? CGLO_PRESSURE : 0)
1211 | (bPres ? CGLO_CONSTRAINT : 0)
1212 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1213 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1216 /* explanation of above:
1217 a) We compute Ekin at the full time step
1218 if 1) we are using the AveVel Ekin, and it's not the
1219 initial step, or 2) if we are using AveEkin, but need the full
1220 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1221 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1222 EkinAveVel because it's needed for the pressure */
1224 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1229 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1230 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1237 /* We need the kinetic energy at minus the half step for determining
1238 * the full step kinetic energy and possibly for T-coupling.*/
1239 /* This may not be quite working correctly yet . . . . */
1240 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1241 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1242 constr, NULL, FALSE, state->box,
1243 top_global, &bSumEkinhOld,
1244 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1249 if (iterate.bIterationActive &&
1250 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1251 state->veta, &vetanew))
1255 bFirstIterate = FALSE;
1258 if (bTrotter && !bInitStep)
1260 copy_mat(shake_vir, state->svir_prev);
1261 copy_mat(force_vir, state->fvir_prev);
1262 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1264 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1265 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1266 enerd->term[F_EKIN] = trace(ekind->ekin);
1269 /* if it's the initial step, we performed this first step just to get the constraint virial */
1270 if (bInitStep && ir->eI == eiVV)
1272 copy_rvecn(cbuf, state->v, 0, state->natoms);
1276 /* MRS -- now done iterating -- compute the conserved quantity */
1279 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1282 last_ekin = enerd->term[F_EKIN];
1284 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1286 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1288 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1291 sum_dhdl(enerd, state->lambda, ir->fepvals);
1295 /* ######## END FIRST UPDATE STEP ############## */
1296 /* ######## If doing VV, we now have v(dt) ###### */
1299 /* perform extended ensemble sampling in lambda - we don't
1300 actually move to the new state before outputting
1301 statistics, but if performing simulated tempering, we
1302 do update the velocities and the tau_t. */
1304 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1305 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1306 copy_df_history(&state_global->dfhist, &state->dfhist);
1309 /* Now we have the energies and forces corresponding to the
1310 * coordinates at time t. We must output all of this before
1313 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1314 ir, state, state_global, top_global, fr,
1315 outf, mdebin, ekind, f, f_global,
1317 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1319 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1320 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1322 /* kludge -- virial is lost with restart for NPT control. Must restart */
1323 if (bStartingFromCpt && bVV)
1325 copy_mat(state->svir_prev, shake_vir);
1326 copy_mat(state->fvir_prev, force_vir);
1329 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1331 /* Check whether everything is still allright */
1332 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1333 #ifdef GMX_THREAD_MPI
1338 /* this is just make gs.sig compatible with the hack
1339 of sending signals around by MPI_Reduce with together with
1341 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1343 gs.sig[eglsSTOPCOND] = 1;
1345 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1347 gs.sig[eglsSTOPCOND] = -1;
1349 /* < 0 means stop at next step, > 0 means stop at next NS step */
1353 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1354 gmx_get_signal_name(),
1355 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1359 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1360 gmx_get_signal_name(),
1361 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1363 handled_stop_condition = (int)gmx_get_stop_condition();
1365 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1366 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1367 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1369 /* Signal to terminate the run */
1370 gs.sig[eglsSTOPCOND] = 1;
1373 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1375 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1378 if (bResetCountersHalfMaxH && MASTER(cr) &&
1379 elapsed_time > max_hours*60.0*60.0*0.495)
1381 gs.sig[eglsRESETCOUNTERS] = 1;
1384 if (ir->nstlist == -1 && !bRerunMD)
1386 /* When bGStatEveryStep=FALSE, global_stat is only called
1387 * when we check the atom displacements, not at NS steps.
1388 * This means that also the bonded interaction count check is not
1389 * performed immediately after NS. Therefore a few MD steps could
1390 * be performed with missing interactions.
1391 * But wrong energies are never written to file,
1392 * since energies are only written after global_stat
1395 if (step >= nlh.step_nscheck)
1397 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1398 nlh.scale_tot, state->x);
1402 /* This is not necessarily true,
1403 * but step_nscheck is determined quite conservatively.
1409 /* In parallel we only have to check for checkpointing in steps
1410 * where we do global communication,
1411 * otherwise the other nodes don't know.
1413 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1416 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1417 gs.set[eglsCHKPT] == 0)
1419 gs.sig[eglsCHKPT] = 1;
1422 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1427 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1429 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1431 gmx_bool bIfRandomize;
1432 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1433 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1434 if (constr && bIfRandomize)
1436 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1437 state, fr->bMolPBC, graph, f,
1438 &top->idef, tmp_vir,
1439 cr, nrnb, wcycle, upd, constr,
1440 TRUE, bCalcVir, vetanew);
1445 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1447 gmx_iterate_init(&iterate, TRUE);
1448 /* for iterations, we save these vectors, as we will be redoing the calculations */
1449 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1452 bFirstIterate = TRUE;
1453 while (bFirstIterate || iterate.bIterationActive)
1455 /* We now restore these vectors to redo the calculation with improved extended variables */
1456 if (iterate.bIterationActive)
1458 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1461 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1462 so scroll down for that logic */
1464 /* ######### START SECOND UPDATE STEP ################# */
1465 /* Box is changed in update() when we do pressure coupling,
1466 * but we should still use the old box for energy corrections and when
1467 * writing it to the energy file, so it matches the trajectory files for
1468 * the same timestep above. Make a copy in a separate array.
1470 copy_mat(state->box, lastbox);
1475 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1477 wallcycle_start(wcycle, ewcUPDATE);
1478 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1481 if (iterate.bIterationActive)
1489 /* we use a new value of scalevir to converge the iterations faster */
1490 scalevir = tracevir/trace(shake_vir);
1492 msmul(shake_vir, scalevir, shake_vir);
1493 m_add(force_vir, shake_vir, total_vir);
1494 clear_mat(shake_vir);
1496 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1497 /* We can only do Berendsen coupling after we have summed
1498 * the kinetic energy or virial. Since the happens
1499 * in global_state after update, we should only do it at
1500 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1505 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1506 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1511 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1513 /* velocity half-step update */
1514 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1515 bUpdateDoLR, fr->f_twin, fcd,
1516 ekind, M, upd, FALSE, etrtVELOCITY2,
1517 cr, nrnb, constr, &top->idef);
1520 /* Above, initialize just copies ekinh into ekin,
1521 * it doesn't copy position (for VV),
1522 * and entire integrator for MD.
1525 if (ir->eI == eiVVAK)
1527 copy_rvecn(state->x, cbuf, 0, state->natoms);
1529 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1531 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1532 bUpdateDoLR, fr->f_twin, fcd,
1533 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1534 wallcycle_stop(wcycle, ewcUPDATE);
1536 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1537 fr->bMolPBC, graph, f,
1538 &top->idef, shake_vir,
1539 cr, nrnb, wcycle, upd, constr,
1540 FALSE, bCalcVir, state->veta);
1542 if (ir->eI == eiVVAK)
1544 /* erase F_EKIN and F_TEMP here? */
1545 /* just compute the kinetic energy at the half step to perform a trotter step */
1546 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1547 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1548 constr, NULL, FALSE, lastbox,
1549 top_global, &bSumEkinhOld,
1550 cglo_flags | CGLO_TEMPERATURE
1552 wallcycle_start(wcycle, ewcUPDATE);
1553 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1554 /* now we know the scaling, we can compute the positions again again */
1555 copy_rvecn(cbuf, state->x, 0, state->natoms);
1557 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1559 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1560 bUpdateDoLR, fr->f_twin, fcd,
1561 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1562 wallcycle_stop(wcycle, ewcUPDATE);
1564 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1565 /* are the small terms in the shake_vir here due
1566 * to numerical errors, or are they important
1567 * physically? I'm thinking they are just errors, but not completely sure.
1568 * For now, will call without actually constraining, constr=NULL*/
1569 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1570 state, fr->bMolPBC, graph, f,
1571 &top->idef, tmp_vir,
1572 cr, nrnb, wcycle, upd, NULL,
1578 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1581 if (fr->bSepDVDL && fplog && do_log)
1583 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1587 /* this factor or 2 correction is necessary
1588 because half of the constraint force is removed
1589 in the vv step, so we have to double it. See
1590 the Redmine issue #1255. It is not yet clear
1591 if the factor of 2 is exact, or just a very
1592 good approximation, and this will be
1593 investigated. The next step is to see if this
1594 can be done adding a dhdl contribution from the
1595 rattle step, but this is somewhat more
1596 complicated with the current code. Will be
1597 investigated, hopefully for 4.6.3. However,
1598 this current solution is much better than
1599 having it completely wrong.
1601 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1605 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1610 /* Need to unshift here */
1611 unshift_self(graph, state->box, state->x);
1616 wallcycle_start(wcycle, ewcVSITECONSTR);
1619 shift_self(graph, state->box, state->x);
1621 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1622 top->idef.iparams, top->idef.il,
1623 fr->ePBC, fr->bMolPBC, cr, state->box);
1627 unshift_self(graph, state->box, state->x);
1629 wallcycle_stop(wcycle, ewcVSITECONSTR);
1632 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1633 /* With Leap-Frog we can skip compute_globals at
1634 * non-communication steps, but we need to calculate
1635 * the kinetic energy one step before communication.
1637 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1639 if (ir->nstlist == -1 && bFirstIterate)
1641 gs.sig[eglsNABNSB] = nlh.nabnsb;
1643 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1644 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1646 bFirstIterate ? &gs : NULL,
1647 (step_rel % gs.nstms == 0) &&
1648 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1650 top_global, &bSumEkinhOld,
1652 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1653 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1654 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1655 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1656 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1657 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1660 if (ir->nstlist == -1 && bFirstIterate)
1662 nlh.nabnsb = gs.set[eglsNABNSB];
1663 gs.set[eglsNABNSB] = 0;
1666 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1667 /* ############# END CALC EKIN AND PRESSURE ################# */
1669 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1670 the virial that should probably be addressed eventually. state->veta has better properies,
1671 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1672 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1674 if (iterate.bIterationActive &&
1675 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1676 trace(shake_vir), &tracevir))
1680 bFirstIterate = FALSE;
1683 if (!bVV || bRerunMD)
1685 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1686 sum_dhdl(enerd, state->lambda, ir->fepvals);
1688 update_box(fplog, step, ir, mdatoms, state, f,
1689 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1691 /* ################# END UPDATE STEP 2 ################# */
1692 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1694 /* The coordinates (x) were unshifted in update */
1697 /* We will not sum ekinh_old,
1698 * so signal that we still have to do it.
1700 bSumEkinhOld = TRUE;
1703 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1705 /* use the directly determined last velocity, not actually the averaged half steps */
1706 if (bTrotter && ir->eI == eiVV)
1708 enerd->term[F_EKIN] = last_ekin;
1710 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1714 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1718 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1720 /* ######### END PREPARING EDR OUTPUT ########### */
1725 gmx_bool do_dr, do_or;
1727 if (fplog && do_log && bDoExpanded)
1729 /* only needed if doing expanded ensemble */
1730 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1731 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1733 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1737 upd_mdebin(mdebin, bDoDHDL, TRUE,
1738 t, mdatoms->tmass, enerd, state,
1739 ir->fepvals, ir->expandedvals, lastbox,
1740 shake_vir, force_vir, total_vir, pres,
1741 ekind, mu_tot, constr);
1745 upd_mdebin_step(mdebin);
1748 do_dr = do_per_step(step, ir->nstdisreout);
1749 do_or = do_per_step(step, ir->nstorireout);
1751 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1753 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1755 if (ir->ePull != epullNO)
1757 pull_print_output(ir->pull, step, t);
1760 if (do_per_step(step, ir->nstlog))
1762 if (fflush(fplog) != 0)
1764 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1770 /* Have to do this part _after_ outputting the logfile and the edr file */
1771 /* Gets written into the state at the beginning of next loop*/
1772 state->fep_state = lamnew;
1774 /* Print the remaining wall clock time for the run */
1775 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1779 fprintf(stderr, "\n");
1781 print_time(stderr, walltime_accounting, step, ir, cr);
1784 /* Ion/water position swapping.
1785 * Not done in last step since trajectory writing happens before this call
1786 * in the MD loop and exchanges would be lost anyway. */
1787 bNeedRepartition = FALSE;
1788 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1789 do_per_step(step, ir->swap->nstswap))
1791 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1792 bRerunMD ? rerun_fr.x : state->x,
1793 bRerunMD ? rerun_fr.box : state->box,
1794 top_global, MASTER(cr) && bVerbose, bRerunMD);
1796 if (bNeedRepartition && DOMAINDECOMP(cr))
1798 dd_collect_state(cr->dd, state, state_global);
1802 /* Replica exchange */
1804 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1805 do_per_step(step, repl_ex_nst))
1807 bExchanged = replica_exchange(fplog, cr, repl_ex,
1808 state_global, enerd,
1812 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1814 dd_partition_system(fplog, step, cr, TRUE, 1,
1815 state_global, top_global, ir,
1816 state, &f, mdatoms, top, fr,
1817 vsite, shellfc, constr,
1818 nrnb, wcycle, FALSE);
1823 bStartingFromCpt = FALSE;
1825 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1826 /* With all integrators, except VV, we need to retain the pressure
1827 * at the current step for coupling at the next step.
1829 if ((state->flags & (1<<estPRES_PREV)) &&
1831 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1833 /* Store the pressure in t_state for pressure coupling
1834 * at the next MD step.
1836 copy_mat(pres, state->pres_prev);
1839 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1841 if ( (membed != NULL) && (!bLastStep) )
1843 rescale_membed(step_rel, membed, state_global->x);
1850 /* read next frame from input trajectory */
1851 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1856 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1860 if (!bRerunMD || !rerun_fr.bStep)
1862 /* increase the MD step number */
1867 cycles = wallcycle_stop(wcycle, ewcSTEP);
1868 if (DOMAINDECOMP(cr) && wcycle)
1870 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1873 if (bPMETuneRunning || bPMETuneTry)
1875 /* PME grid + cut-off optimization with GPUs or PME nodes */
1877 /* Count the total cycles over the last steps */
1878 cycles_pmes += cycles;
1880 /* We can only switch cut-off at NS steps */
1881 if (step % ir->nstlist == 0)
1883 /* PME grid + cut-off optimization with GPUs or PME nodes */
1886 if (DDMASTER(cr->dd))
1888 /* PME node load is too high, start tuning */
1889 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1891 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1893 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1895 bPMETuneTry = FALSE;
1898 if (bPMETuneRunning)
1900 /* init_step might not be a multiple of nstlist,
1901 * but the first cycle is always skipped anyhow.
1904 pme_load_balance(pme_loadbal, cr,
1905 (bVerbose && MASTER(cr)) ? stderr : NULL,
1907 ir, state, cycles_pmes,
1908 fr->ic, fr->nbv, &fr->pmedata,
1911 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1912 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1913 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1914 fr->rlist = fr->ic->rlist;
1915 fr->rlistlong = fr->ic->rlistlong;
1916 fr->rcoulomb = fr->ic->rcoulomb;
1917 fr->rvdw = fr->ic->rvdw;
1919 if (ir->eDispCorr != edispcNO)
1921 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1928 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1929 gs.set[eglsRESETCOUNTERS] != 0)
1931 /* Reset all the counters related to performance over the run */
1932 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1933 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1934 wcycle_set_reset_counters(wcycle, -1);
1935 if (!(cr->duty & DUTY_PME))
1937 /* Tell our PME node to reset its counters */
1938 gmx_pme_send_resetcounters(cr, step);
1940 /* Correct max_hours for the elapsed time */
1941 max_hours -= elapsed_time/(60.0*60.0);
1942 bResetCountersHalfMaxH = FALSE;
1943 gs.set[eglsRESETCOUNTERS] = 0;
1946 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1947 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1950 /* End of main MD loop */
1953 /* Stop measuring walltime */
1954 walltime_accounting_end(walltime_accounting);
1956 if (bRerunMD && MASTER(cr))
1961 if (!(cr->duty & DUTY_PME))
1963 /* Tell the PME only node to finish */
1964 gmx_pme_send_finish(cr);
1969 if (ir->nstcalcenergy > 0 && !bRerunMD)
1971 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1972 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1979 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1981 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)));
1982 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1985 if (pme_loadbal != NULL)
1987 pme_loadbal_done(pme_loadbal, cr, fplog,
1988 fr->nbv != NULL && fr->nbv->bUseGPU);
1991 if (shellfc && fplog)
1993 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1994 (nconverged*100.0)/step_rel);
1995 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1999 if (repl_ex_nst > 0 && MASTER(cr))
2001 print_replica_exchange_statistics(fplog, repl_ex);
2004 /* IMD cleanup, if bIMD is TRUE. */
2005 IMD_finalize(ir->bIMD, ir->imd);
2007 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);