2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
42 #include "gromacs/utility/smalloc.h"
54 #include "md_support.h"
55 #include "md_logging.h"
69 #include "domdec_network.h"
70 #include "gromacs/gmxlib/topsort.h"
74 #include "gromacs/gmxpreprocess/compute_io.h"
75 #include "checkpoint.h"
76 #include "mtop_util.h"
77 #include "sighandler.h"
79 #include "gromacs/utility/cstringutil.h"
80 #include "pme_loadbal.h"
83 #include "types/nlistheuristics.h"
84 #include "types/iteratedconstraints.h"
85 #include "nbnxn_cuda_data_mgmt.h"
87 #include "gromacs/utility/gmxmpi.h"
88 #include "gromacs/fileio/confio.h"
89 #include "gromacs/fileio/trajectory_writing.h"
90 #include "gromacs/fileio/trnio.h"
91 #include "gromacs/fileio/trxio.h"
92 #include "gromacs/fileio/xtcio.h"
93 #include "gromacs/timing/wallcycle.h"
94 #include "gromacs/timing/walltime_accounting.h"
95 #include "gromacs/pulling/pull.h"
96 #include "gromacs/swap/swapcoords.h"
97 #include "gromacs/imd/imd.h"
100 #include "corewrap.h"
103 static void reset_all_counters(FILE *fplog, t_commrec *cr,
105 gmx_int64_t *step_rel, t_inputrec *ir,
106 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
107 gmx_walltime_accounting_t walltime_accounting,
108 nbnxn_cuda_ptr_t cu_nbv)
110 char sbuf[STEPSTRSIZE];
112 /* Reset all the counters related to performance over the run */
113 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
114 gmx_step_str(step, sbuf));
118 nbnxn_cuda_reset_timings(cu_nbv);
121 wallcycle_stop(wcycle, ewcRUN);
122 wallcycle_reset_all(wcycle);
123 if (DOMAINDECOMP(cr))
125 reset_dd_statistics_counters(cr->dd);
128 ir->init_step += *step_rel;
129 ir->nsteps -= *step_rel;
131 wallcycle_start(wcycle, ewcRUN);
132 walltime_accounting_start(walltime_accounting);
133 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
136 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
137 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
139 gmx_vsite_t *vsite, gmx_constr_t constr,
140 int stepout, t_inputrec *ir,
141 gmx_mtop_t *top_global,
143 t_state *state_global,
145 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
146 gmx_edsam_t ed, t_forcerec *fr,
147 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
148 real cpt_period, real max_hours,
149 const char gmx_unused *deviceOptions,
152 gmx_walltime_accounting_t walltime_accounting)
154 gmx_mdoutf_t outf = NULL;
155 gmx_int64_t step, step_rel;
157 double t, t0, lam0[efptNR];
158 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
159 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
160 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
161 bBornRadii, bStartingFromCpt;
162 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
163 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
164 bForceUpdate = FALSE, bCPT;
165 gmx_bool bMasterState;
166 int force_flags, cglo_flags;
167 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
172 t_state *bufstate = NULL;
173 matrix *scale_tot, pcoupl_mu, M, ebox;
176 gmx_repl_ex_t repl_ex = NULL;
179 t_mdebin *mdebin = NULL;
180 t_state *state = NULL;
181 rvec *f_global = NULL;
182 gmx_enerdata_t *enerd;
184 gmx_global_stat_t gstat;
185 gmx_update_t upd = NULL;
186 t_graph *graph = NULL;
188 gmx_groups_t *groups;
189 gmx_ekindata_t *ekind, *ekind_save;
190 gmx_shellfc_t shellfc;
191 int count, nconverged = 0;
194 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged, bNeedRepartition;
196 gmx_bool bResetCountersHalfMaxH = FALSE;
197 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
198 gmx_bool bUpdateDoLR;
202 real veta_save, scalevir, tracevir;
208 real saved_conserved_quantity = 0;
213 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
214 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
215 gmx_iterate_t iterate;
216 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
217 simulation stops. If equal to zero, don't
218 communicate any more between multisims.*/
219 /* PME load balancing data for GPU kernels */
220 pme_load_balancing_t pme_loadbal = NULL;
222 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
225 gmx_bool bIMDstep = FALSE;
228 /* Temporary addition for FAHCORE checkpointing */
232 /* Check for special mdrun options */
233 bRerunMD = (Flags & MD_RERUN);
234 bAppend = (Flags & MD_APPENDFILES);
235 if (Flags & MD_RESETCOUNTERSHALFWAY)
239 /* Signal to reset the counters half the simulation steps. */
240 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
242 /* Signal to reset the counters halfway the simulation time. */
243 bResetCountersHalfMaxH = (max_hours > 0);
246 /* md-vv uses averaged full step velocities for T-control
247 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
248 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
250 if (bVV) /* to store the initial velocities while computing virial */
252 snew(cbuf, top_global->natoms);
254 /* all the iteratative cases - only if there are constraints */
255 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
256 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
257 false in this step. The correct value, true or false,
258 is set at each step, as it depends on the frequency of temperature
259 and pressure control.*/
260 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
264 /* Since we don't know if the frames read are related in any way,
265 * rebuild the neighborlist at every step.
268 ir->nstcalcenergy = 1;
272 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
274 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
275 bGStatEveryStep = (nstglobalcomm == 1);
277 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
280 "To reduce the energy communication with nstlist = -1\n"
281 "the neighbor list validity should not be checked at every step,\n"
282 "this means that exact integration is not guaranteed.\n"
283 "The neighbor list validity is checked after:\n"
284 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
285 "In most cases this will result in exact integration.\n"
286 "This reduces the energy communication by a factor of 2 to 3.\n"
287 "If you want less energy communication, set nstlist > 3.\n\n");
292 ir->nstxout_compressed = 0;
294 groups = &top_global->groups;
297 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
298 &(state_global->fep_state), lam0,
299 nrnb, top_global, &upd,
300 nfile, fnm, &outf, &mdebin,
301 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
303 clear_mat(total_vir);
305 /* Energy terms and groups */
307 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
309 if (DOMAINDECOMP(cr))
315 snew(f, top_global->natoms);
318 /* Kinetic energy data */
320 init_ekindata(fplog, top_global, &(ir->opts), ekind);
321 /* needed for iteration of constraints */
323 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
324 /* Copy the cos acceleration to the groups struct */
325 ekind->cosacc.cos_accel = ir->cos_accel;
327 gstat = global_stat_init(ir);
330 /* Check for polarizable models and flexible constraints */
331 shellfc = init_shell_flexcon(fplog, fr->cutoff_scheme == ecutsVERLET,
332 top_global, n_flexible_constraints(constr),
333 (ir->bContinuation ||
334 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
335 NULL : state_global->x);
339 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
340 set_deform_reference_box(upd,
341 deform_init_init_step_tpx,
342 deform_init_box_tpx);
343 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
347 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
348 if ((io > 2000) && MASTER(cr))
351 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
356 if (DOMAINDECOMP(cr))
358 top = dd_init_local_top(top_global);
361 dd_init_local_state(cr->dd, state_global, state);
363 if (DDMASTER(cr->dd) && ir->nstfout)
365 snew(f_global, state_global->natoms);
370 top = gmx_mtop_generate_local_top(top_global, ir);
372 forcerec_set_excl_load(fr, top);
374 state = serial_init_local_state(state_global);
377 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
381 set_vsite_top(vsite, top, mdatoms, cr);
384 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
386 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
391 make_local_shells(cr, mdatoms, shellfc);
394 setup_bonded_threading(fr, &top->idef);
397 /* Set up interactive MD (IMD) */
398 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
399 nfile, fnm, oenv, imdport, Flags);
401 if (DOMAINDECOMP(cr))
403 /* Distribute the charge groups over the nodes from the master node */
404 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
405 state_global, top_global, ir,
406 state, &f, mdatoms, top, fr,
407 vsite, shellfc, constr,
408 nrnb, wcycle, FALSE);
412 update_mdatoms(mdatoms, state->lambda[efptMASS]);
414 if (opt2bSet("-cpi", nfile, fnm))
416 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
420 bStateFromCP = FALSE;
425 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
432 /* Update mdebin with energy history if appending to output files */
433 if (Flags & MD_APPENDFILES)
435 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
439 /* We might have read an energy history from checkpoint,
440 * free the allocated memory and reset the counts.
442 done_energyhistory(&state_global->enerhist);
443 init_energyhistory(&state_global->enerhist);
446 /* Set the initial energy history in state by updating once */
447 update_energyhistory(&state_global->enerhist, mdebin);
450 /* Initialize constraints */
451 if (constr && !DOMAINDECOMP(cr))
453 set_constraints(constr, top, ir, mdatoms, cr);
458 /* We need to be sure replica exchange can only occur
459 * when the energies are current */
460 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
461 "repl_ex_nst", &repl_ex_nst);
462 /* This check needs to happen before inter-simulation
463 * signals are initialized, too */
465 if (repl_ex_nst > 0 && MASTER(cr))
467 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
468 repl_ex_nst, repl_ex_nex, repl_ex_seed);
471 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
472 * PME tuning is not supported with PME only for LJ and not for Coulomb.
474 if ((Flags & MD_TUNEPME) &&
475 EEL_PME(fr->eeltype) &&
476 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
479 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
481 if (cr->duty & DUTY_PME)
483 /* Start tuning right away, as we can't measure the load */
484 bPMETuneRunning = TRUE;
488 /* Separate PME nodes, we can measure the PP/PME load balance */
493 if (!ir->bContinuation && !bRerunMD)
495 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
497 /* Set the velocities of frozen particles to zero */
498 for (i = 0; i < mdatoms->homenr; i++)
500 for (m = 0; m < DIM; m++)
502 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
512 /* Constrain the initial coordinates and velocities */
513 do_constrain_first(fplog, constr, ir, mdatoms, state,
518 /* Construct the virtual sites for the initial configuration */
519 construct_vsites(vsite, state->x, ir->delta_t, NULL,
520 top->idef.iparams, top->idef.il,
521 fr->ePBC, fr->bMolPBC, cr, state->box);
527 /* set free energy calculation frequency as the minimum
528 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
529 nstfep = ir->fepvals->nstdhdl;
532 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
536 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
539 /* I'm assuming we need global communication the first time! MRS */
540 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
541 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
542 | (bVV ? CGLO_PRESSURE : 0)
543 | (bVV ? CGLO_CONSTRAINT : 0)
544 | (bRerunMD ? CGLO_RERUNMD : 0)
545 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
547 bSumEkinhOld = FALSE;
548 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
549 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
550 constr, NULL, FALSE, state->box,
551 top_global, &bSumEkinhOld, cglo_flags);
552 if (ir->eI == eiVVAK)
554 /* a second call to get the half step temperature initialized as well */
555 /* we do the same call as above, but turn the pressure off -- internally to
556 compute_globals, this is recognized as a velocity verlet half-step
557 kinetic energy calculation. This minimized excess variables, but
558 perhaps loses some logic?*/
560 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
561 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
562 constr, NULL, FALSE, state->box,
563 top_global, &bSumEkinhOld,
564 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
567 /* Calculate the initial half step temperature, and save the ekinh_old */
568 if (!(Flags & MD_STARTFROMCPT))
570 for (i = 0; (i < ir->opts.ngtc); i++)
572 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
577 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
578 and there is no previous step */
581 /* if using an iterative algorithm, we need to create a working directory for the state. */
584 bufstate = init_bufstate(state);
587 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
588 temperature control */
589 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
593 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
596 "RMS relative constraint deviation after constraining: %.2e\n",
597 constr_rmsd(constr, FALSE));
599 if (EI_STATE_VELOCITY(ir->eI))
601 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
605 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
606 " input trajectory '%s'\n\n",
607 *(top_global->name), opt2fn("-rerun", nfile, fnm));
610 fprintf(stderr, "Calculated time to finish depends on nsteps from "
611 "run input file,\nwhich may not correspond to the time "
612 "needed to process input trajectory.\n\n");
618 fprintf(stderr, "starting mdrun '%s'\n",
619 *(top_global->name));
622 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
626 sprintf(tbuf, "%s", "infinite");
628 if (ir->init_step > 0)
630 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
631 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
632 gmx_step_str(ir->init_step, sbuf2),
633 ir->init_step*ir->delta_t);
637 fprintf(stderr, "%s steps, %s ps.\n",
638 gmx_step_str(ir->nsteps, sbuf), tbuf);
641 fprintf(fplog, "\n");
644 walltime_accounting_start(walltime_accounting);
645 wallcycle_start(wcycle, ewcRUN);
646 print_start(fplog, cr, walltime_accounting, "mdrun");
648 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
650 chkpt_ret = fcCheckPointParallel( cr->nodeid,
654 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
659 /***********************************************************
663 ************************************************************/
665 /* if rerunMD then read coordinates and velocities from input trajectory */
668 if (getenv("GMX_FORCE_UPDATE"))
676 bNotLastFrame = read_first_frame(oenv, &status,
677 opt2fn("-rerun", nfile, fnm),
678 &rerun_fr, TRX_NEED_X | TRX_READ_V);
679 if (rerun_fr.natoms != top_global->natoms)
682 "Number of atoms in trajectory (%d) does not match the "
683 "run input file (%d)\n",
684 rerun_fr.natoms, top_global->natoms);
686 if (ir->ePBC != epbcNONE)
690 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);
692 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
694 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
701 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
704 if (ir->ePBC != epbcNONE)
706 /* Set the shift vectors.
707 * Necessary here when have a static box different from the tpr box.
709 calc_shifts(rerun_fr.box, fr->shift_vec);
713 /* loop over MD steps or if rerunMD to end of input trajectory */
715 /* Skip the first Nose-Hoover integration when we get the state from tpx */
716 bStateFromTPX = !bStateFromCP;
717 bInitStep = bFirstStep && (bStateFromTPX || bVV);
718 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
720 bSumEkinhOld = FALSE;
722 bNeedRepartition = FALSE;
724 init_global_signals(&gs, cr, ir, repl_ex_nst);
726 step = ir->init_step;
729 if (ir->nstlist == -1)
731 init_nlistheuristics(&nlh, bGStatEveryStep, step);
734 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
736 /* check how many steps are left in other sims */
737 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
741 /* and stop now if we should */
742 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
743 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
744 while (!bLastStep || (bRerunMD && bNotLastFrame))
747 wallcycle_start(wcycle, ewcSTEP);
753 step = rerun_fr.step;
754 step_rel = step - ir->init_step;
767 bLastStep = (step_rel == ir->nsteps);
768 t = t0 + step*ir->delta_t;
771 if (ir->efep != efepNO || ir->bSimTemp)
773 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
774 requiring different logic. */
776 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
777 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
778 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
779 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
780 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
785 update_annealing_target_temp(&(ir->opts), t);
790 if (!DOMAINDECOMP(cr) || MASTER(cr))
792 for (i = 0; i < state_global->natoms; i++)
794 copy_rvec(rerun_fr.x[i], state_global->x[i]);
798 for (i = 0; i < state_global->natoms; i++)
800 copy_rvec(rerun_fr.v[i], state_global->v[i]);
805 for (i = 0; i < state_global->natoms; i++)
807 clear_rvec(state_global->v[i]);
811 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
812 " Ekin, temperature and pressure are incorrect,\n"
813 " the virial will be incorrect when constraints are present.\n"
815 bRerunWarnNoV = FALSE;
819 copy_mat(rerun_fr.box, state_global->box);
820 copy_mat(state_global->box, state->box);
822 if (vsite && (Flags & MD_RERUN_VSITE))
824 if (DOMAINDECOMP(cr))
826 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
830 /* Following is necessary because the graph may get out of sync
831 * with the coordinates if we only have every N'th coordinate set
833 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
834 shift_self(graph, state->box, state->x);
836 construct_vsites(vsite, state->x, ir->delta_t, state->v,
837 top->idef.iparams, top->idef.il,
838 fr->ePBC, fr->bMolPBC, cr, state->box);
841 unshift_self(graph, state->box, state->x);
846 /* Stop Center of Mass motion */
847 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
851 /* for rerun MD always do Neighbour Searching */
852 bNS = (bFirstStep || ir->nstlist != 0);
857 /* Determine whether or not to do Neighbour Searching and LR */
858 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
860 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
861 (ir->nstlist == -1 && nlh.nabnsb > 0));
863 if (bNS && ir->nstlist == -1)
865 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
869 /* check whether we should stop because another simulation has
873 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
874 (multisim_nsteps != ir->nsteps) )
881 "Stopping simulation %d because another one has finished\n",
885 gs.sig[eglsCHKPT] = 1;
890 /* < 0 means stop at next step, > 0 means stop at next NS step */
891 if ( (gs.set[eglsSTOPCOND] < 0) ||
892 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
897 /* Determine whether or not to update the Born radii if doing GB */
898 bBornRadii = bFirstStep;
899 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
904 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
905 do_verbose = bVerbose &&
906 (step % stepout == 0 || bFirstStep || bLastStep);
908 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
916 bMasterState = FALSE;
917 /* Correct the new box if it is too skewed */
918 if (DYNAMIC_BOX(*ir))
920 if (correct_box(fplog, step, state->box, graph))
925 if (DOMAINDECOMP(cr) && bMasterState)
927 dd_collect_state(cr->dd, state, state_global);
931 if (DOMAINDECOMP(cr))
933 /* Repartition the domain decomposition */
934 wallcycle_start(wcycle, ewcDOMDEC);
935 dd_partition_system(fplog, step, cr,
936 bMasterState, nstglobalcomm,
937 state_global, top_global, ir,
938 state, &f, mdatoms, top, fr,
939 vsite, shellfc, constr,
941 do_verbose && !bPMETuneRunning);
942 wallcycle_stop(wcycle, ewcDOMDEC);
943 /* If using an iterative integrator, reallocate space to match the decomposition */
947 if (MASTER(cr) && do_log)
949 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
952 if (ir->efep != efepNO)
954 update_mdatoms(mdatoms, state->lambda[efptMASS]);
957 if ((bRerunMD && rerun_fr.bV) || bExchanged)
960 /* We need the kinetic energy at minus the half step for determining
961 * the full step kinetic energy and possibly for T-coupling.*/
962 /* This may not be quite working correctly yet . . . . */
963 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
964 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
965 constr, NULL, FALSE, state->box,
966 top_global, &bSumEkinhOld,
967 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
969 clear_mat(force_vir);
971 /* We write a checkpoint at this MD step when:
972 * either at an NS step when we signalled through gs,
973 * or at the last step (but not when we do not want confout),
974 * but never at the first step or with rerun.
976 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
977 (bLastStep && (Flags & MD_CONFOUT))) &&
978 step > ir->init_step && !bRerunMD);
981 gs.set[eglsCHKPT] = 0;
984 /* Determine the energy and pressure:
985 * at nstcalcenergy steps and at energy output steps (set below).
987 if (EI_VV(ir->eI) && (!bInitStep))
989 /* for vv, the first half of the integration actually corresponds
990 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
991 but the virial needs to be calculated on both the current step and the 'next' step. Future
992 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
994 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
995 bCalcVir = bCalcEner ||
996 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1000 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1001 bCalcVir = bCalcEner ||
1002 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1005 /* Do we need global communication ? */
1006 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1007 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1008 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1010 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1012 if (do_ene || do_log)
1019 /* these CGLO_ options remain the same throughout the iteration */
1020 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1021 (bGStat ? CGLO_GSTAT : 0)
1024 force_flags = (GMX_FORCE_STATECHANGED |
1025 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1026 GMX_FORCE_ALLFORCES |
1028 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1029 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1030 (bDoFEP ? GMX_FORCE_DHDL : 0)
1035 if (do_per_step(step, ir->nstcalclr))
1037 force_flags |= GMX_FORCE_DO_LR;
1043 /* Now is the time to relax the shells */
1044 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1045 ir, bNS, force_flags,
1048 state, f, force_vir, mdatoms,
1049 nrnb, wcycle, graph, groups,
1050 shellfc, fr, bBornRadii, t, mu_tot,
1052 mdoutf_get_fp_field(outf));
1062 /* The coordinates (x) are shifted (to get whole molecules)
1064 * This is parallellized as well, and does communication too.
1065 * Check comments in sim_util.c
1067 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1068 state->box, state->x, &state->hist,
1069 f, force_vir, mdatoms, enerd, fcd,
1070 state->lambda, graph,
1071 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1072 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1075 if (bVV && !bStartingFromCpt && !bRerunMD)
1076 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1078 if (ir->eI == eiVV && bInitStep)
1080 /* if using velocity verlet with full time step Ekin,
1081 * take the first half step only to compute the
1082 * virial for the first step. From there,
1083 * revert back to the initial coordinates
1084 * so that the input is actually the initial step.
1086 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1090 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1091 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1094 /* If we are using twin-range interactions where the long-range component
1095 * is only evaluated every nstcalclr>1 steps, we should do a special update
1096 * step to combine the long-range forces on these steps.
1097 * For nstcalclr=1 this is not done, since the forces would have been added
1098 * directly to the short-range forces already.
1100 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1102 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1103 f, bUpdateDoLR, fr->f_twin, fcd,
1104 ekind, M, upd, bInitStep, etrtVELOCITY1,
1105 cr, nrnb, constr, &top->idef);
1107 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1109 gmx_iterate_init(&iterate, TRUE);
1111 /* for iterations, we save these vectors, as we will be self-consistently iterating
1114 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1116 /* save the state */
1117 if (iterate.bIterationActive)
1119 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1122 bFirstIterate = TRUE;
1123 while (bFirstIterate || iterate.bIterationActive)
1125 if (iterate.bIterationActive)
1127 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1128 if (bFirstIterate && bTrotter)
1130 /* The first time through, we need a decent first estimate
1131 of veta(t+dt) to compute the constraints. Do
1132 this by computing the box volume part of the
1133 trotter integration at this time. Nothing else
1134 should be changed by this routine here. If
1135 !(first time), we start with the previous value
1138 veta_save = state->veta;
1139 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1140 vetanew = state->veta;
1141 state->veta = veta_save;
1146 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1148 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1149 state, fr->bMolPBC, graph, f,
1150 &top->idef, shake_vir,
1151 cr, nrnb, wcycle, upd, constr,
1152 TRUE, bCalcVir, vetanew);
1156 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1162 /* Need to unshift here if a do_force has been
1163 called in the previous step */
1164 unshift_self(graph, state->box, state->x);
1167 /* if VV, compute the pressure and constraints */
1168 /* For VV2, we strictly only need this if using pressure
1169 * control, but we really would like to have accurate pressures
1171 * Think about ways around this in the future?
1172 * For now, keep this choice in comments.
1174 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1175 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1177 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1178 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1180 bSumEkinhOld = TRUE;
1182 /* for vv, the first half of the integration actually corresponds to the previous step.
1183 So we need information from the last step in the first half of the integration */
1184 if (bGStat || do_per_step(step-1, nstglobalcomm))
1186 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1187 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1188 constr, NULL, FALSE, state->box,
1189 top_global, &bSumEkinhOld,
1192 | (bTemp ? CGLO_TEMPERATURE : 0)
1193 | (bPres ? CGLO_PRESSURE : 0)
1194 | (bPres ? CGLO_CONSTRAINT : 0)
1195 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1196 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1199 /* explanation of above:
1200 a) We compute Ekin at the full time step
1201 if 1) we are using the AveVel Ekin, and it's not the
1202 initial step, or 2) if we are using AveEkin, but need the full
1203 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1204 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1205 EkinAveVel because it's needed for the pressure */
1207 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1212 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1213 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1220 /* We need the kinetic energy at minus the half step for determining
1221 * the full step kinetic energy and possibly for T-coupling.*/
1222 /* This may not be quite working correctly yet . . . . */
1223 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1224 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1225 constr, NULL, FALSE, state->box,
1226 top_global, &bSumEkinhOld,
1227 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1232 if (iterate.bIterationActive &&
1233 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1234 state->veta, &vetanew))
1238 bFirstIterate = FALSE;
1241 if (bTrotter && !bInitStep)
1243 copy_mat(shake_vir, state->svir_prev);
1244 copy_mat(force_vir, state->fvir_prev);
1245 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1247 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1248 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1249 enerd->term[F_EKIN] = trace(ekind->ekin);
1252 /* if it's the initial step, we performed this first step just to get the constraint virial */
1253 if (bInitStep && ir->eI == eiVV)
1255 copy_rvecn(cbuf, state->v, 0, state->natoms);
1259 /* MRS -- now done iterating -- compute the conserved quantity */
1262 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1265 last_ekin = enerd->term[F_EKIN];
1267 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1269 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1271 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1274 sum_dhdl(enerd, state->lambda, ir->fepvals);
1278 /* ######## END FIRST UPDATE STEP ############## */
1279 /* ######## If doing VV, we now have v(dt) ###### */
1282 /* perform extended ensemble sampling in lambda - we don't
1283 actually move to the new state before outputting
1284 statistics, but if performing simulated tempering, we
1285 do update the velocities and the tau_t. */
1287 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1288 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1289 copy_df_history(&state_global->dfhist, &state->dfhist);
1292 /* Now we have the energies and forces corresponding to the
1293 * coordinates at time t. We must output all of this before
1296 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1297 ir, state, state_global, top_global, fr,
1298 outf, mdebin, ekind, f, f_global,
1300 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1302 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1303 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1305 /* kludge -- virial is lost with restart for NPT control. Must restart */
1306 if (bStartingFromCpt && bVV)
1308 copy_mat(state->svir_prev, shake_vir);
1309 copy_mat(state->fvir_prev, force_vir);
1312 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1314 /* Check whether everything is still allright */
1315 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1316 #ifdef GMX_THREAD_MPI
1321 /* this is just make gs.sig compatible with the hack
1322 of sending signals around by MPI_Reduce with together with
1324 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1326 gs.sig[eglsSTOPCOND] = 1;
1328 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1330 gs.sig[eglsSTOPCOND] = -1;
1332 /* < 0 means stop at next step, > 0 means stop at next NS step */
1336 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1337 gmx_get_signal_name(),
1338 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1342 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1343 gmx_get_signal_name(),
1344 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1346 handled_stop_condition = (int)gmx_get_stop_condition();
1348 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1349 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1350 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1352 /* Signal to terminate the run */
1353 gs.sig[eglsSTOPCOND] = 1;
1356 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1358 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1361 if (bResetCountersHalfMaxH && MASTER(cr) &&
1362 elapsed_time > max_hours*60.0*60.0*0.495)
1364 gs.sig[eglsRESETCOUNTERS] = 1;
1367 if (ir->nstlist == -1 && !bRerunMD)
1369 /* When bGStatEveryStep=FALSE, global_stat is only called
1370 * when we check the atom displacements, not at NS steps.
1371 * This means that also the bonded interaction count check is not
1372 * performed immediately after NS. Therefore a few MD steps could
1373 * be performed with missing interactions.
1374 * But wrong energies are never written to file,
1375 * since energies are only written after global_stat
1378 if (step >= nlh.step_nscheck)
1380 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1381 nlh.scale_tot, state->x);
1385 /* This is not necessarily true,
1386 * but step_nscheck is determined quite conservatively.
1392 /* In parallel we only have to check for checkpointing in steps
1393 * where we do global communication,
1394 * otherwise the other nodes don't know.
1396 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1399 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1400 gs.set[eglsCHKPT] == 0)
1402 gs.sig[eglsCHKPT] = 1;
1405 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1410 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1412 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1414 gmx_bool bIfRandomize;
1415 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1416 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1417 if (constr && bIfRandomize)
1419 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1420 state, fr->bMolPBC, graph, f,
1421 &top->idef, tmp_vir,
1422 cr, nrnb, wcycle, upd, constr,
1423 TRUE, bCalcVir, vetanew);
1428 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1430 gmx_iterate_init(&iterate, TRUE);
1431 /* for iterations, we save these vectors, as we will be redoing the calculations */
1432 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1435 bFirstIterate = TRUE;
1436 while (bFirstIterate || iterate.bIterationActive)
1438 /* We now restore these vectors to redo the calculation with improved extended variables */
1439 if (iterate.bIterationActive)
1441 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1444 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1445 so scroll down for that logic */
1447 /* ######### START SECOND UPDATE STEP ################# */
1448 /* Box is changed in update() when we do pressure coupling,
1449 * but we should still use the old box for energy corrections and when
1450 * writing it to the energy file, so it matches the trajectory files for
1451 * the same timestep above. Make a copy in a separate array.
1453 copy_mat(state->box, lastbox);
1458 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1460 wallcycle_start(wcycle, ewcUPDATE);
1461 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1464 if (iterate.bIterationActive)
1472 /* we use a new value of scalevir to converge the iterations faster */
1473 scalevir = tracevir/trace(shake_vir);
1475 msmul(shake_vir, scalevir, shake_vir);
1476 m_add(force_vir, shake_vir, total_vir);
1477 clear_mat(shake_vir);
1479 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1480 /* We can only do Berendsen coupling after we have summed
1481 * the kinetic energy or virial. Since the happens
1482 * in global_state after update, we should only do it at
1483 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1488 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1489 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1494 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1496 /* velocity half-step update */
1497 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1498 bUpdateDoLR, fr->f_twin, fcd,
1499 ekind, M, upd, FALSE, etrtVELOCITY2,
1500 cr, nrnb, constr, &top->idef);
1503 /* Above, initialize just copies ekinh into ekin,
1504 * it doesn't copy position (for VV),
1505 * and entire integrator for MD.
1508 if (ir->eI == eiVVAK)
1510 copy_rvecn(state->x, cbuf, 0, state->natoms);
1512 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1514 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1515 bUpdateDoLR, fr->f_twin, fcd,
1516 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1517 wallcycle_stop(wcycle, ewcUPDATE);
1519 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1520 fr->bMolPBC, graph, f,
1521 &top->idef, shake_vir,
1522 cr, nrnb, wcycle, upd, constr,
1523 FALSE, bCalcVir, state->veta);
1525 if (ir->eI == eiVVAK)
1527 /* erase F_EKIN and F_TEMP here? */
1528 /* just compute the kinetic energy at the half step to perform a trotter step */
1529 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1530 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1531 constr, NULL, FALSE, lastbox,
1532 top_global, &bSumEkinhOld,
1533 cglo_flags | CGLO_TEMPERATURE
1535 wallcycle_start(wcycle, ewcUPDATE);
1536 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1537 /* now we know the scaling, we can compute the positions again again */
1538 copy_rvecn(cbuf, state->x, 0, state->natoms);
1540 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1542 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1543 bUpdateDoLR, fr->f_twin, fcd,
1544 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1545 wallcycle_stop(wcycle, ewcUPDATE);
1547 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1548 /* are the small terms in the shake_vir here due
1549 * to numerical errors, or are they important
1550 * physically? I'm thinking they are just errors, but not completely sure.
1551 * For now, will call without actually constraining, constr=NULL*/
1552 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1553 state, fr->bMolPBC, graph, f,
1554 &top->idef, tmp_vir,
1555 cr, nrnb, wcycle, upd, NULL,
1561 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1564 if (fr->bSepDVDL && fplog && do_log)
1566 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1570 /* this factor or 2 correction is necessary
1571 because half of the constraint force is removed
1572 in the vv step, so we have to double it. See
1573 the Redmine issue #1255. It is not yet clear
1574 if the factor of 2 is exact, or just a very
1575 good approximation, and this will be
1576 investigated. The next step is to see if this
1577 can be done adding a dhdl contribution from the
1578 rattle step, but this is somewhat more
1579 complicated with the current code. Will be
1580 investigated, hopefully for 4.6.3. However,
1581 this current solution is much better than
1582 having it completely wrong.
1584 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1588 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1593 /* Need to unshift here */
1594 unshift_self(graph, state->box, state->x);
1599 wallcycle_start(wcycle, ewcVSITECONSTR);
1602 shift_self(graph, state->box, state->x);
1604 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1605 top->idef.iparams, top->idef.il,
1606 fr->ePBC, fr->bMolPBC, cr, state->box);
1610 unshift_self(graph, state->box, state->x);
1612 wallcycle_stop(wcycle, ewcVSITECONSTR);
1615 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1616 /* With Leap-Frog we can skip compute_globals at
1617 * non-communication steps, but we need to calculate
1618 * the kinetic energy one step before communication.
1620 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1622 if (ir->nstlist == -1 && bFirstIterate)
1624 gs.sig[eglsNABNSB] = nlh.nabnsb;
1626 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1627 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1629 bFirstIterate ? &gs : NULL,
1630 (step_rel % gs.nstms == 0) &&
1631 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1633 top_global, &bSumEkinhOld,
1635 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1636 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1637 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1638 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1639 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1640 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1643 if (ir->nstlist == -1 && bFirstIterate)
1645 nlh.nabnsb = gs.set[eglsNABNSB];
1646 gs.set[eglsNABNSB] = 0;
1649 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1650 /* ############# END CALC EKIN AND PRESSURE ################# */
1652 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1653 the virial that should probably be addressed eventually. state->veta has better properies,
1654 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1655 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1657 if (iterate.bIterationActive &&
1658 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1659 trace(shake_vir), &tracevir))
1663 bFirstIterate = FALSE;
1666 if (!bVV || bRerunMD)
1668 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1669 sum_dhdl(enerd, state->lambda, ir->fepvals);
1671 update_box(fplog, step, ir, mdatoms, state, f,
1672 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1674 /* ################# END UPDATE STEP 2 ################# */
1675 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1677 /* The coordinates (x) were unshifted in update */
1680 /* We will not sum ekinh_old,
1681 * so signal that we still have to do it.
1683 bSumEkinhOld = TRUE;
1686 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1688 /* use the directly determined last velocity, not actually the averaged half steps */
1689 if (bTrotter && ir->eI == eiVV)
1691 enerd->term[F_EKIN] = last_ekin;
1693 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1697 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1701 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1703 /* ######### END PREPARING EDR OUTPUT ########### */
1708 gmx_bool do_dr, do_or;
1710 if (fplog && do_log && bDoExpanded)
1712 /* only needed if doing expanded ensemble */
1713 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1714 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1716 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1720 upd_mdebin(mdebin, bDoDHDL, TRUE,
1721 t, mdatoms->tmass, enerd, state,
1722 ir->fepvals, ir->expandedvals, lastbox,
1723 shake_vir, force_vir, total_vir, pres,
1724 ekind, mu_tot, constr);
1728 upd_mdebin_step(mdebin);
1731 do_dr = do_per_step(step, ir->nstdisreout);
1732 do_or = do_per_step(step, ir->nstorireout);
1734 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1736 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1738 if (ir->ePull != epullNO)
1740 pull_print_output(ir->pull, step, t);
1743 if (do_per_step(step, ir->nstlog))
1745 if (fflush(fplog) != 0)
1747 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1753 /* Have to do this part _after_ outputting the logfile and the edr file */
1754 /* Gets written into the state at the beginning of next loop*/
1755 state->fep_state = lamnew;
1757 /* Print the remaining wall clock time for the run */
1758 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1762 fprintf(stderr, "\n");
1764 print_time(stderr, walltime_accounting, step, ir, cr);
1767 /* Ion/water position swapping.
1768 * Not done in last step since trajectory writing happens before this call
1769 * in the MD loop and exchanges would be lost anyway. */
1770 bNeedRepartition = FALSE;
1771 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1772 do_per_step(step, ir->swap->nstswap))
1774 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1775 bRerunMD ? rerun_fr.x : state->x,
1776 bRerunMD ? rerun_fr.box : state->box,
1777 top_global, MASTER(cr) && bVerbose, bRerunMD);
1779 if (bNeedRepartition && DOMAINDECOMP(cr))
1781 dd_collect_state(cr->dd, state, state_global);
1785 /* Replica exchange */
1787 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1788 do_per_step(step, repl_ex_nst))
1790 bExchanged = replica_exchange(fplog, cr, repl_ex,
1791 state_global, enerd,
1795 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1797 dd_partition_system(fplog, step, cr, TRUE, 1,
1798 state_global, top_global, ir,
1799 state, &f, mdatoms, top, fr,
1800 vsite, shellfc, constr,
1801 nrnb, wcycle, FALSE);
1806 bStartingFromCpt = FALSE;
1808 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1809 /* With all integrators, except VV, we need to retain the pressure
1810 * at the current step for coupling at the next step.
1812 if ((state->flags & (1<<estPRES_PREV)) &&
1814 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1816 /* Store the pressure in t_state for pressure coupling
1817 * at the next MD step.
1819 copy_mat(pres, state->pres_prev);
1822 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1824 if ( (membed != NULL) && (!bLastStep) )
1826 rescale_membed(step_rel, membed, state_global->x);
1833 /* read next frame from input trajectory */
1834 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1839 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1843 if (!bRerunMD || !rerun_fr.bStep)
1845 /* increase the MD step number */
1850 cycles = wallcycle_stop(wcycle, ewcSTEP);
1851 if (DOMAINDECOMP(cr) && wcycle)
1853 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1856 if (bPMETuneRunning || bPMETuneTry)
1858 /* PME grid + cut-off optimization with GPUs or PME nodes */
1860 /* Count the total cycles over the last steps */
1861 cycles_pmes += cycles;
1863 /* We can only switch cut-off at NS steps */
1864 if (step % ir->nstlist == 0)
1866 /* PME grid + cut-off optimization with GPUs or PME nodes */
1869 if (DDMASTER(cr->dd))
1871 /* PME node load is too high, start tuning */
1872 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1874 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1876 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1878 bPMETuneTry = FALSE;
1881 if (bPMETuneRunning)
1883 /* init_step might not be a multiple of nstlist,
1884 * but the first cycle is always skipped anyhow.
1887 pme_load_balance(pme_loadbal, cr,
1888 (bVerbose && MASTER(cr)) ? stderr : NULL,
1890 ir, state, cycles_pmes,
1891 fr->ic, fr->nbv, &fr->pmedata,
1894 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1895 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1896 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1897 fr->rlist = fr->ic->rlist;
1898 fr->rlistlong = fr->ic->rlistlong;
1899 fr->rcoulomb = fr->ic->rcoulomb;
1900 fr->rvdw = fr->ic->rvdw;
1902 if (ir->eDispCorr != edispcNO)
1904 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1911 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1912 gs.set[eglsRESETCOUNTERS] != 0)
1914 /* Reset all the counters related to performance over the run */
1915 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1916 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1917 wcycle_set_reset_counters(wcycle, -1);
1918 if (!(cr->duty & DUTY_PME))
1920 /* Tell our PME node to reset its counters */
1921 gmx_pme_send_resetcounters(cr, step);
1923 /* Correct max_hours for the elapsed time */
1924 max_hours -= elapsed_time/(60.0*60.0);
1925 bResetCountersHalfMaxH = FALSE;
1926 gs.set[eglsRESETCOUNTERS] = 0;
1929 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1930 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1933 /* End of main MD loop */
1936 /* Stop measuring walltime */
1937 walltime_accounting_end(walltime_accounting);
1939 if (bRerunMD && MASTER(cr))
1944 if (!(cr->duty & DUTY_PME))
1946 /* Tell the PME only node to finish */
1947 gmx_pme_send_finish(cr);
1952 if (ir->nstcalcenergy > 0 && !bRerunMD)
1954 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1955 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1962 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1964 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)));
1965 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1968 if (pme_loadbal != NULL)
1970 pme_loadbal_done(pme_loadbal, cr, fplog,
1971 fr->nbv != NULL && fr->nbv->bUseGPU);
1974 if (shellfc && fplog)
1976 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1977 (nconverged*100.0)/step_rel);
1978 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1982 if (repl_ex_nst > 0 && MASTER(cr))
1984 print_replica_exchange_statistics(fplog, repl_ex);
1987 /* IMD cleanup, if bIMD is TRUE. */
1988 IMD_finalize(ir->bIMD, ir->imd);
1990 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);