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/math/vec.h"
51 #include "md_support.h"
52 #include "md_logging.h"
64 #include "domdec_network.h"
68 #include "gromacs/gmxpreprocess/compute_io.h"
69 #include "checkpoint.h"
70 #include "gromacs/topology/mtop_util.h"
71 #include "sighandler.h"
73 #include "gromacs/utility/cstringutil.h"
74 #include "pme_loadbal.h"
77 #include "types/nlistheuristics.h"
78 #include "types/iteratedconstraints.h"
79 #include "nbnxn_cuda_data_mgmt.h"
81 #include "gromacs/fileio/confio.h"
82 #include "gromacs/fileio/mdoutf.h"
83 #include "gromacs/fileio/trajectory_writing.h"
84 #include "gromacs/fileio/trnio.h"
85 #include "gromacs/fileio/trxio.h"
86 #include "gromacs/fileio/xtcio.h"
87 #include "gromacs/imd/imd.h"
88 #include "gromacs/pbcutil/mshift.h"
89 #include "gromacs/pbcutil/pbc.h"
90 #include "gromacs/pulling/pull.h"
91 #include "gromacs/swap/swapcoords.h"
92 #include "gromacs/timing/wallcycle.h"
93 #include "gromacs/timing/walltime_accounting.h"
94 #include "gromacs/utility/gmxmpi.h"
95 #include "gromacs/utility/smalloc.h"
101 static void reset_all_counters(FILE *fplog, t_commrec *cr,
103 gmx_int64_t *step_rel, t_inputrec *ir,
104 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
105 gmx_walltime_accounting_t walltime_accounting,
106 nbnxn_cuda_ptr_t cu_nbv)
108 char sbuf[STEPSTRSIZE];
110 /* Reset all the counters related to performance over the run */
111 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
112 gmx_step_str(step, sbuf));
116 nbnxn_cuda_reset_timings(cu_nbv);
119 wallcycle_stop(wcycle, ewcRUN);
120 wallcycle_reset_all(wcycle);
121 if (DOMAINDECOMP(cr))
123 reset_dd_statistics_counters(cr->dd);
126 ir->init_step += *step_rel;
127 ir->nsteps -= *step_rel;
129 wallcycle_start(wcycle, ewcRUN);
130 walltime_accounting_start(walltime_accounting);
131 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
134 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
135 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
137 gmx_vsite_t *vsite, gmx_constr_t constr,
138 int stepout, t_inputrec *ir,
139 gmx_mtop_t *top_global,
141 t_state *state_global,
143 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
144 gmx_edsam_t ed, t_forcerec *fr,
145 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
146 real cpt_period, real max_hours,
147 const char gmx_unused *deviceOptions,
150 gmx_walltime_accounting_t walltime_accounting)
152 gmx_mdoutf_t outf = NULL;
153 gmx_int64_t step, step_rel;
155 double t, t0, lam0[efptNR];
156 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
157 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
158 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
159 bBornRadii, bStartingFromCpt;
160 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
161 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
162 bForceUpdate = FALSE, bCPT;
163 gmx_bool bMasterState;
164 int force_flags, cglo_flags;
165 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
170 t_state *bufstate = NULL;
174 gmx_repl_ex_t repl_ex = NULL;
177 t_mdebin *mdebin = NULL;
178 t_state *state = NULL;
179 rvec *f_global = NULL;
180 gmx_enerdata_t *enerd;
182 gmx_global_stat_t gstat;
183 gmx_update_t upd = NULL;
184 t_graph *graph = NULL;
186 gmx_groups_t *groups;
187 gmx_ekindata_t *ekind, *ekind_save;
188 gmx_shellfc_t shellfc;
189 int count, nconverged = 0;
191 gmx_bool bConverged = TRUE, bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
192 gmx_bool bResetCountersHalfMaxH = FALSE;
193 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
194 gmx_bool bUpdateDoLR;
198 real veta_save, scalevir, tracevir;
204 real saved_conserved_quantity = 0;
208 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
209 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
210 gmx_iterate_t iterate;
211 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
212 simulation stops. If equal to zero, don't
213 communicate any more between multisims.*/
214 /* PME load balancing data for GPU kernels */
215 pme_load_balancing_t pme_loadbal = NULL;
217 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
220 gmx_bool bIMDstep = FALSE;
223 /* Temporary addition for FAHCORE checkpointing */
227 /* Check for special mdrun options */
228 bRerunMD = (Flags & MD_RERUN);
229 if (Flags & MD_RESETCOUNTERSHALFWAY)
233 /* Signal to reset the counters half the simulation steps. */
234 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
236 /* Signal to reset the counters halfway the simulation time. */
237 bResetCountersHalfMaxH = (max_hours > 0);
240 /* md-vv uses averaged full step velocities for T-control
241 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
242 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
244 if (bVV) /* to store the initial velocities while computing virial */
246 snew(cbuf, top_global->natoms);
248 /* all the iteratative cases - only if there are constraints */
249 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
250 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
251 false in this step. The correct value, true or false,
252 is set at each step, as it depends on the frequency of temperature
253 and pressure control.*/
254 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
258 /* Since we don't know if the frames read are related in any way,
259 * rebuild the neighborlist at every step.
262 ir->nstcalcenergy = 1;
266 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
268 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
269 bGStatEveryStep = (nstglobalcomm == 1);
271 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
274 "To reduce the energy communication with nstlist = -1\n"
275 "the neighbor list validity should not be checked at every step,\n"
276 "this means that exact integration is not guaranteed.\n"
277 "The neighbor list validity is checked after:\n"
278 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
279 "In most cases this will result in exact integration.\n"
280 "This reduces the energy communication by a factor of 2 to 3.\n"
281 "If you want less energy communication, set nstlist > 3.\n\n");
286 ir->nstxout_compressed = 0;
288 groups = &top_global->groups;
291 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
292 &(state_global->fep_state), lam0,
293 nrnb, top_global, &upd,
294 nfile, fnm, &outf, &mdebin,
295 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
297 clear_mat(total_vir);
299 /* Energy terms and groups */
301 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
303 if (DOMAINDECOMP(cr))
309 snew(f, top_global->natoms);
312 /* Kinetic energy data */
314 init_ekindata(fplog, top_global, &(ir->opts), ekind);
315 /* needed for iteration of constraints */
317 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
318 /* Copy the cos acceleration to the groups struct */
319 ekind->cosacc.cos_accel = ir->cos_accel;
321 gstat = global_stat_init(ir);
324 /* Check for polarizable models and flexible constraints */
325 shellfc = init_shell_flexcon(fplog,
326 top_global, n_flexible_constraints(constr),
327 (ir->bContinuation ||
328 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
329 NULL : state_global->x);
330 if (shellfc && ir->nstcalcenergy != 1)
332 gmx_fatal(FARGS, "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is not supported in combinations with shell particles.\nPlease make a new tpr file.", ir->nstcalcenergy);
334 if (shellfc && DOMAINDECOMP(cr))
336 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
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);
469 if (repl_ex_nst > 0 && MASTER(cr))
471 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
472 repl_ex_nst, repl_ex_nex, repl_ex_seed);
475 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
476 * PME tuning is not supported with PME only for LJ and not for Coulomb.
478 if ((Flags & MD_TUNEPME) &&
479 EEL_PME(fr->eeltype) &&
480 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
483 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
485 if (cr->duty & DUTY_PME)
487 /* Start tuning right away, as we can't measure the load */
488 bPMETuneRunning = TRUE;
492 /* Separate PME nodes, we can measure the PP/PME load balance */
497 if (!ir->bContinuation && !bRerunMD)
499 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
501 /* Set the velocities of frozen particles to zero */
502 for (i = 0; i < mdatoms->homenr; i++)
504 for (m = 0; m < DIM; m++)
506 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
516 /* Constrain the initial coordinates and velocities */
517 do_constrain_first(fplog, constr, ir, mdatoms, state,
522 /* Construct the virtual sites for the initial configuration */
523 construct_vsites(vsite, state->x, ir->delta_t, NULL,
524 top->idef.iparams, top->idef.il,
525 fr->ePBC, fr->bMolPBC, cr, state->box);
531 /* set free energy calculation frequency as the minimum
532 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
533 nstfep = ir->fepvals->nstdhdl;
536 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
540 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
543 /* I'm assuming we need global communication the first time! MRS */
544 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
545 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
546 | (bVV ? CGLO_PRESSURE : 0)
547 | (bVV ? CGLO_CONSTRAINT : 0)
548 | (bRerunMD ? CGLO_RERUNMD : 0)
549 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
551 bSumEkinhOld = FALSE;
552 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
553 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
554 constr, NULL, FALSE, state->box,
555 top_global, &bSumEkinhOld, cglo_flags);
556 if (ir->eI == eiVVAK)
558 /* a second call to get the half step temperature initialized as well */
559 /* we do the same call as above, but turn the pressure off -- internally to
560 compute_globals, this is recognized as a velocity verlet half-step
561 kinetic energy calculation. This minimized excess variables, but
562 perhaps loses some logic?*/
564 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
565 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
566 constr, NULL, FALSE, state->box,
567 top_global, &bSumEkinhOld,
568 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
571 /* Calculate the initial half step temperature, and save the ekinh_old */
572 if (!(Flags & MD_STARTFROMCPT))
574 for (i = 0; (i < ir->opts.ngtc); i++)
576 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
581 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
582 and there is no previous step */
585 /* if using an iterative algorithm, we need to create a working directory for the state. */
588 bufstate = init_bufstate(state);
591 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
592 temperature control */
593 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
597 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
600 "RMS relative constraint deviation after constraining: %.2e\n",
601 constr_rmsd(constr, FALSE));
603 if (EI_STATE_VELOCITY(ir->eI))
605 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
609 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
610 " input trajectory '%s'\n\n",
611 *(top_global->name), opt2fn("-rerun", nfile, fnm));
614 fprintf(stderr, "Calculated time to finish depends on nsteps from "
615 "run input file,\nwhich may not correspond to the time "
616 "needed to process input trajectory.\n\n");
622 fprintf(stderr, "starting mdrun '%s'\n",
623 *(top_global->name));
626 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
630 sprintf(tbuf, "%s", "infinite");
632 if (ir->init_step > 0)
634 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
635 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
636 gmx_step_str(ir->init_step, sbuf2),
637 ir->init_step*ir->delta_t);
641 fprintf(stderr, "%s steps, %s ps.\n",
642 gmx_step_str(ir->nsteps, sbuf), tbuf);
645 fprintf(fplog, "\n");
648 walltime_accounting_start(walltime_accounting);
649 wallcycle_start(wcycle, ewcRUN);
650 print_start(fplog, cr, walltime_accounting, "mdrun");
652 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
654 chkpt_ret = fcCheckPointParallel( cr->nodeid,
658 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
663 /***********************************************************
667 ************************************************************/
669 /* if rerunMD then read coordinates and velocities from input trajectory */
672 if (getenv("GMX_FORCE_UPDATE"))
680 bNotLastFrame = read_first_frame(oenv, &status,
681 opt2fn("-rerun", nfile, fnm),
682 &rerun_fr, TRX_NEED_X | TRX_READ_V);
683 if (rerun_fr.natoms != top_global->natoms)
686 "Number of atoms in trajectory (%d) does not match the "
687 "run input file (%d)\n",
688 rerun_fr.natoms, top_global->natoms);
690 if (ir->ePBC != epbcNONE)
694 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);
696 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
698 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
705 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
708 if (ir->ePBC != epbcNONE)
710 /* Set the shift vectors.
711 * Necessary here when have a static box different from the tpr box.
713 calc_shifts(rerun_fr.box, fr->shift_vec);
717 /* loop over MD steps or if rerunMD to end of input trajectory */
719 /* Skip the first Nose-Hoover integration when we get the state from tpx */
720 bStateFromTPX = !bStateFromCP;
721 bInitStep = bFirstStep && (bStateFromTPX || bVV);
722 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
723 bSumEkinhOld = FALSE;
725 bNeedRepartition = FALSE;
727 init_global_signals(&gs, cr, ir, repl_ex_nst);
729 step = ir->init_step;
732 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));
783 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
784 do_per_step(step, repl_ex_nst));
788 update_annealing_target_temp(&(ir->opts), t);
793 if (!DOMAINDECOMP(cr) || MASTER(cr))
795 for (i = 0; i < state_global->natoms; i++)
797 copy_rvec(rerun_fr.x[i], state_global->x[i]);
801 for (i = 0; i < state_global->natoms; i++)
803 copy_rvec(rerun_fr.v[i], state_global->v[i]);
808 for (i = 0; i < state_global->natoms; i++)
810 clear_rvec(state_global->v[i]);
814 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
815 " Ekin, temperature and pressure are incorrect,\n"
816 " the virial will be incorrect when constraints are present.\n"
818 bRerunWarnNoV = FALSE;
822 copy_mat(rerun_fr.box, state_global->box);
823 copy_mat(state_global->box, state->box);
825 if (vsite && (Flags & MD_RERUN_VSITE))
827 if (DOMAINDECOMP(cr))
829 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
833 /* Following is necessary because the graph may get out of sync
834 * with the coordinates if we only have every N'th coordinate set
836 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
837 shift_self(graph, state->box, state->x);
839 construct_vsites(vsite, state->x, ir->delta_t, state->v,
840 top->idef.iparams, top->idef.il,
841 fr->ePBC, fr->bMolPBC, cr, state->box);
844 unshift_self(graph, state->box, state->x);
849 /* Stop Center of Mass motion */
850 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
854 /* for rerun MD always do Neighbour Searching */
855 bNS = (bFirstStep || ir->nstlist != 0);
860 /* Determine whether or not to do Neighbour Searching and LR */
861 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
863 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
864 (ir->nstlist == -1 && nlh.nabnsb > 0));
866 if (bNS && ir->nstlist == -1)
868 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
872 /* check whether we should stop because another simulation has
876 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
877 (multisim_nsteps != ir->nsteps) )
884 "Stopping simulation %d because another one has finished\n",
888 gs.sig[eglsCHKPT] = 1;
893 /* < 0 means stop at next step, > 0 means stop at next NS step */
894 if ( (gs.set[eglsSTOPCOND] < 0) ||
895 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
900 /* Determine whether or not to update the Born radii if doing GB */
901 bBornRadii = bFirstStep;
902 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
907 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
908 do_verbose = bVerbose &&
909 (step % stepout == 0 || bFirstStep || bLastStep);
911 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
919 bMasterState = FALSE;
920 /* Correct the new box if it is too skewed */
921 if (DYNAMIC_BOX(*ir))
923 if (correct_box(fplog, step, state->box, graph))
928 if (DOMAINDECOMP(cr) && bMasterState)
930 dd_collect_state(cr->dd, state, state_global);
934 if (DOMAINDECOMP(cr))
936 /* Repartition the domain decomposition */
937 wallcycle_start(wcycle, ewcDOMDEC);
938 dd_partition_system(fplog, step, cr,
939 bMasterState, nstglobalcomm,
940 state_global, top_global, ir,
941 state, &f, mdatoms, top, fr,
942 vsite, shellfc, constr,
944 do_verbose && !bPMETuneRunning);
945 wallcycle_stop(wcycle, ewcDOMDEC);
946 /* If using an iterative integrator, reallocate space to match the decomposition */
950 if (MASTER(cr) && do_log)
952 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
955 if (ir->efep != efepNO)
957 update_mdatoms(mdatoms, state->lambda[efptMASS]);
960 if ((bRerunMD && rerun_fr.bV) || bExchanged)
963 /* We need the kinetic energy at minus the half step for determining
964 * the full step kinetic energy and possibly for T-coupling.*/
965 /* This may not be quite working correctly yet . . . . */
966 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
967 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
968 constr, NULL, FALSE, state->box,
969 top_global, &bSumEkinhOld,
970 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
972 clear_mat(force_vir);
974 /* We write a checkpoint at this MD step when:
975 * either at an NS step when we signalled through gs,
976 * or at the last step (but not when we do not want confout),
977 * but never at the first step or with rerun.
979 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
980 (bLastStep && (Flags & MD_CONFOUT))) &&
981 step > ir->init_step && !bRerunMD);
984 gs.set[eglsCHKPT] = 0;
987 /* Determine the energy and pressure:
988 * at nstcalcenergy steps and at energy output steps (set below).
990 if (EI_VV(ir->eI) && (!bInitStep))
992 /* for vv, the first half of the integration actually corresponds
993 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
994 but the virial needs to be calculated on both the current step and the 'next' step. Future
995 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
997 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
998 bCalcVir = bCalcEner ||
999 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1003 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1004 bCalcVir = bCalcEner ||
1005 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1008 /* Do we need global communication ? */
1009 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1010 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1011 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1013 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1015 if (do_ene || do_log || bDoReplEx)
1022 /* these CGLO_ options remain the same throughout the iteration */
1023 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1024 (bGStat ? CGLO_GSTAT : 0)
1027 force_flags = (GMX_FORCE_STATECHANGED |
1028 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1029 GMX_FORCE_ALLFORCES |
1031 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1032 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1033 (bDoFEP ? GMX_FORCE_DHDL : 0)
1038 if (do_per_step(step, ir->nstcalclr))
1040 force_flags |= GMX_FORCE_DO_LR;
1046 /* Now is the time to relax the shells */
1047 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1048 ir, bNS, force_flags,
1051 state, f, force_vir, mdatoms,
1052 nrnb, wcycle, graph, groups,
1053 shellfc, fr, bBornRadii, t, mu_tot,
1055 mdoutf_get_fp_field(outf));
1065 /* The coordinates (x) are shifted (to get whole molecules)
1067 * This is parallellized as well, and does communication too.
1068 * Check comments in sim_util.c
1070 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1071 state->box, state->x, &state->hist,
1072 f, force_vir, mdatoms, enerd, fcd,
1073 state->lambda, graph,
1074 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1075 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1078 if (bVV && !bStartingFromCpt && !bRerunMD)
1079 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1081 if (ir->eI == eiVV && bInitStep)
1083 /* if using velocity verlet with full time step Ekin,
1084 * take the first half step only to compute the
1085 * virial for the first step. From there,
1086 * revert back to the initial coordinates
1087 * so that the input is actually the initial step.
1089 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1093 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1094 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1097 /* If we are using twin-range interactions where the long-range component
1098 * is only evaluated every nstcalclr>1 steps, we should do a special update
1099 * step to combine the long-range forces on these steps.
1100 * For nstcalclr=1 this is not done, since the forces would have been added
1101 * directly to the short-range forces already.
1103 * TODO Remove various aspects of VV+twin-range in master
1104 * branch, because VV integrators did not ever support
1105 * twin-range multiple time stepping with constraints.
1107 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1109 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1110 f, bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1111 ekind, M, upd, bInitStep, etrtVELOCITY1,
1112 cr, nrnb, constr, &top->idef);
1114 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1116 gmx_iterate_init(&iterate, TRUE);
1118 /* for iterations, we save these vectors, as we will be self-consistently iterating
1121 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1123 /* save the state */
1124 if (iterate.bIterationActive)
1126 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1129 bFirstIterate = TRUE;
1130 while (bFirstIterate || iterate.bIterationActive)
1132 if (iterate.bIterationActive)
1134 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1135 if (bFirstIterate && bTrotter)
1137 /* The first time through, we need a decent first estimate
1138 of veta(t+dt) to compute the constraints. Do
1139 this by computing the box volume part of the
1140 trotter integration at this time. Nothing else
1141 should be changed by this routine here. If
1142 !(first time), we start with the previous value
1145 veta_save = state->veta;
1146 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1147 vetanew = state->veta;
1148 state->veta = veta_save;
1152 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1154 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1155 state, fr->bMolPBC, graph, f,
1156 &top->idef, shake_vir,
1157 cr, nrnb, wcycle, upd, constr,
1158 TRUE, bCalcVir, vetanew);
1160 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1162 /* Correct the virial for multiple time stepping */
1163 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1168 /* Need to unshift here if a do_force has been
1169 called in the previous step */
1170 unshift_self(graph, state->box, state->x);
1173 /* if VV, compute the pressure and constraints */
1174 /* For VV2, we strictly only need this if using pressure
1175 * control, but we really would like to have accurate pressures
1177 * Think about ways around this in the future?
1178 * For now, keep this choice in comments.
1180 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1181 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1183 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1184 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1186 bSumEkinhOld = TRUE;
1188 /* for vv, the first half of the integration actually corresponds to the previous step.
1189 So we need information from the last step in the first half of the integration */
1190 if (bGStat || do_per_step(step-1, nstglobalcomm))
1192 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1193 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1194 constr, NULL, FALSE, state->box,
1195 top_global, &bSumEkinhOld,
1198 | (bTemp ? CGLO_TEMPERATURE : 0)
1199 | (bPres ? CGLO_PRESSURE : 0)
1200 | (bPres ? CGLO_CONSTRAINT : 0)
1201 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1202 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1205 /* explanation of above:
1206 a) We compute Ekin at the full time step
1207 if 1) we are using the AveVel Ekin, and it's not the
1208 initial step, or 2) if we are using AveEkin, but need the full
1209 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1210 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1211 EkinAveVel because it's needed for the pressure */
1213 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1218 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1219 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1226 /* We need the kinetic energy at minus the half step for determining
1227 * the full step kinetic energy and possibly for T-coupling.*/
1228 /* This may not be quite working correctly yet . . . . */
1229 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1230 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1231 constr, NULL, FALSE, state->box,
1232 top_global, &bSumEkinhOld,
1233 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1238 if (iterate.bIterationActive &&
1239 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1240 state->veta, &vetanew))
1244 bFirstIterate = FALSE;
1247 if (bTrotter && !bInitStep)
1249 copy_mat(shake_vir, state->svir_prev);
1250 copy_mat(force_vir, state->fvir_prev);
1251 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1253 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1254 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1255 enerd->term[F_EKIN] = trace(ekind->ekin);
1258 /* if it's the initial step, we performed this first step just to get the constraint virial */
1259 if (bInitStep && ir->eI == eiVV)
1261 copy_rvecn(cbuf, state->v, 0, state->natoms);
1265 /* MRS -- now done iterating -- compute the conserved quantity */
1268 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1271 last_ekin = enerd->term[F_EKIN];
1273 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1275 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1277 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1280 sum_dhdl(enerd, state->lambda, ir->fepvals);
1284 /* ######## END FIRST UPDATE STEP ############## */
1285 /* ######## If doing VV, we now have v(dt) ###### */
1288 /* perform extended ensemble sampling in lambda - we don't
1289 actually move to the new state before outputting
1290 statistics, but if performing simulated tempering, we
1291 do update the velocities and the tau_t. */
1293 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1294 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1295 copy_df_history(&state_global->dfhist, &state->dfhist);
1298 /* Now we have the energies and forces corresponding to the
1299 * coordinates at time t. We must output all of this before
1302 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1303 ir, state, state_global, top_global, fr,
1304 outf, mdebin, ekind, f, f_global,
1306 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1308 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1309 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1311 /* kludge -- virial is lost with restart for NPT control. Must restart */
1312 if (bStartingFromCpt && bVV)
1314 copy_mat(state->svir_prev, shake_vir);
1315 copy_mat(state->fvir_prev, force_vir);
1318 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1320 /* Check whether everything is still allright */
1321 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1322 #ifdef GMX_THREAD_MPI
1327 /* this is just make gs.sig compatible with the hack
1328 of sending signals around by MPI_Reduce with together with
1330 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1332 gs.sig[eglsSTOPCOND] = 1;
1334 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1336 gs.sig[eglsSTOPCOND] = -1;
1338 /* < 0 means stop at next step, > 0 means stop at next NS step */
1342 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1343 gmx_get_signal_name(),
1344 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1348 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1349 gmx_get_signal_name(),
1350 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1352 handled_stop_condition = (int)gmx_get_stop_condition();
1354 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1355 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1356 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1358 /* Signal to terminate the run */
1359 gs.sig[eglsSTOPCOND] = 1;
1362 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1364 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1367 if (bResetCountersHalfMaxH && MASTER(cr) &&
1368 elapsed_time > max_hours*60.0*60.0*0.495)
1370 gs.sig[eglsRESETCOUNTERS] = 1;
1373 if (ir->nstlist == -1 && !bRerunMD)
1375 /* When bGStatEveryStep=FALSE, global_stat is only called
1376 * when we check the atom displacements, not at NS steps.
1377 * This means that also the bonded interaction count check is not
1378 * performed immediately after NS. Therefore a few MD steps could
1379 * be performed with missing interactions.
1380 * But wrong energies are never written to file,
1381 * since energies are only written after global_stat
1384 if (step >= nlh.step_nscheck)
1386 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1387 nlh.scale_tot, state->x);
1391 /* This is not necessarily true,
1392 * but step_nscheck is determined quite conservatively.
1398 /* In parallel we only have to check for checkpointing in steps
1399 * where we do global communication,
1400 * otherwise the other nodes don't know.
1402 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1405 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1406 gs.set[eglsCHKPT] == 0)
1408 gs.sig[eglsCHKPT] = 1;
1411 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1416 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1418 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1420 gmx_bool bIfRandomize;
1421 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1422 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1423 if (constr && bIfRandomize)
1425 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1426 state, fr->bMolPBC, graph, f,
1427 &top->idef, tmp_vir,
1428 cr, nrnb, wcycle, upd, constr,
1429 TRUE, bCalcVir, vetanew);
1434 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1436 gmx_iterate_init(&iterate, TRUE);
1437 /* for iterations, we save these vectors, as we will be redoing the calculations */
1438 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1441 bFirstIterate = TRUE;
1442 while (bFirstIterate || iterate.bIterationActive)
1444 /* We now restore these vectors to redo the calculation with improved extended variables */
1445 if (iterate.bIterationActive)
1447 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1450 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1451 so scroll down for that logic */
1453 /* ######### START SECOND UPDATE STEP ################# */
1454 /* Box is changed in update() when we do pressure coupling,
1455 * but we should still use the old box for energy corrections and when
1456 * writing it to the energy file, so it matches the trajectory files for
1457 * the same timestep above. Make a copy in a separate array.
1459 copy_mat(state->box, lastbox);
1463 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1465 wallcycle_start(wcycle, ewcUPDATE);
1466 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1469 if (iterate.bIterationActive)
1477 /* we use a new value of scalevir to converge the iterations faster */
1478 scalevir = tracevir/trace(shake_vir);
1480 msmul(shake_vir, scalevir, shake_vir);
1481 m_add(force_vir, shake_vir, total_vir);
1482 clear_mat(shake_vir);
1484 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1485 /* We can only do Berendsen coupling after we have summed
1486 * the kinetic energy or virial. Since the happens
1487 * in global_state after update, we should only do it at
1488 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1493 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1494 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1499 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1501 /* velocity half-step update */
1502 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1503 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1504 ekind, M, upd, FALSE, etrtVELOCITY2,
1505 cr, nrnb, constr, &top->idef);
1508 /* Above, initialize just copies ekinh into ekin,
1509 * it doesn't copy position (for VV),
1510 * and entire integrator for MD.
1513 if (ir->eI == eiVVAK)
1515 copy_rvecn(state->x, cbuf, 0, state->natoms);
1517 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1519 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1520 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1521 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1522 wallcycle_stop(wcycle, ewcUPDATE);
1524 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1525 fr->bMolPBC, graph, f,
1526 &top->idef, shake_vir,
1527 cr, nrnb, wcycle, upd, constr,
1528 FALSE, bCalcVir, state->veta);
1530 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1532 /* Correct the virial for multiple time stepping */
1533 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1536 if (ir->eI == eiVVAK)
1538 /* erase F_EKIN and F_TEMP here? */
1539 /* just compute the kinetic energy at the half step to perform a trotter step */
1540 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1541 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1542 constr, NULL, FALSE, lastbox,
1543 top_global, &bSumEkinhOld,
1544 cglo_flags | CGLO_TEMPERATURE
1546 wallcycle_start(wcycle, ewcUPDATE);
1547 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1548 /* now we know the scaling, we can compute the positions again again */
1549 copy_rvecn(cbuf, state->x, 0, state->natoms);
1551 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1553 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1554 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1555 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1556 wallcycle_stop(wcycle, ewcUPDATE);
1558 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1559 /* are the small terms in the shake_vir here due
1560 * to numerical errors, or are they important
1561 * physically? I'm thinking they are just errors, but not completely sure.
1562 * For now, will call without actually constraining, constr=NULL*/
1563 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1564 state, fr->bMolPBC, graph, f,
1565 &top->idef, tmp_vir,
1566 cr, nrnb, wcycle, upd, NULL,
1571 if (fr->bSepDVDL && fplog && do_log)
1573 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1577 /* this factor or 2 correction is necessary
1578 because half of the constraint force is removed
1579 in the vv step, so we have to double it. See
1580 the Redmine issue #1255. It is not yet clear
1581 if the factor of 2 is exact, or just a very
1582 good approximation, and this will be
1583 investigated. The next step is to see if this
1584 can be done adding a dhdl contribution from the
1585 rattle step, but this is somewhat more
1586 complicated with the current code. Will be
1587 investigated, hopefully for 4.6.3. However,
1588 this current solution is much better than
1589 having it completely wrong.
1591 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1595 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1600 /* Need to unshift here */
1601 unshift_self(graph, state->box, state->x);
1606 wallcycle_start(wcycle, ewcVSITECONSTR);
1609 shift_self(graph, state->box, state->x);
1611 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1612 top->idef.iparams, top->idef.il,
1613 fr->ePBC, fr->bMolPBC, cr, state->box);
1617 unshift_self(graph, state->box, state->x);
1619 wallcycle_stop(wcycle, ewcVSITECONSTR);
1622 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1623 /* With Leap-Frog we can skip compute_globals at
1624 * non-communication steps, but we need to calculate
1625 * the kinetic energy one step before communication.
1627 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1629 if (ir->nstlist == -1 && bFirstIterate)
1631 gs.sig[eglsNABNSB] = nlh.nabnsb;
1633 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1634 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1636 bFirstIterate ? &gs : NULL,
1637 (step_rel % gs.nstms == 0) &&
1638 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1640 top_global, &bSumEkinhOld,
1642 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1643 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1644 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1645 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1646 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1647 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1650 if (ir->nstlist == -1 && bFirstIterate)
1652 nlh.nabnsb = gs.set[eglsNABNSB];
1653 gs.set[eglsNABNSB] = 0;
1656 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1657 /* ############# END CALC EKIN AND PRESSURE ################# */
1659 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1660 the virial that should probably be addressed eventually. state->veta has better properies,
1661 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1662 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1664 if (iterate.bIterationActive &&
1665 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1666 trace(shake_vir), &tracevir))
1670 bFirstIterate = FALSE;
1673 if (!bVV || bRerunMD)
1675 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1676 sum_dhdl(enerd, state->lambda, ir->fepvals);
1678 update_box(fplog, step, ir, mdatoms, state, f,
1679 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1681 /* ################# END UPDATE STEP 2 ################# */
1682 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1684 /* The coordinates (x) were unshifted in update */
1687 /* We will not sum ekinh_old,
1688 * so signal that we still have to do it.
1690 bSumEkinhOld = TRUE;
1693 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1695 /* use the directly determined last velocity, not actually the averaged half steps */
1696 if (bTrotter && ir->eI == eiVV)
1698 enerd->term[F_EKIN] = last_ekin;
1700 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1704 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1708 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1710 /* ######### END PREPARING EDR OUTPUT ########### */
1715 gmx_bool do_dr, do_or;
1717 if (fplog && do_log && bDoExpanded)
1719 /* only needed if doing expanded ensemble */
1720 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1721 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1723 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1727 upd_mdebin(mdebin, bDoDHDL, TRUE,
1728 t, mdatoms->tmass, enerd, state,
1729 ir->fepvals, ir->expandedvals, lastbox,
1730 shake_vir, force_vir, total_vir, pres,
1731 ekind, mu_tot, constr);
1735 upd_mdebin_step(mdebin);
1738 do_dr = do_per_step(step, ir->nstdisreout);
1739 do_or = do_per_step(step, ir->nstorireout);
1741 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1743 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1745 if (ir->ePull != epullNO)
1747 pull_print_output(ir->pull, step, t);
1750 if (do_per_step(step, ir->nstlog))
1752 if (fflush(fplog) != 0)
1754 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1760 /* Have to do this part _after_ outputting the logfile and the edr file */
1761 /* Gets written into the state at the beginning of next loop*/
1762 state->fep_state = lamnew;
1764 /* Print the remaining wall clock time for the run */
1765 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1769 fprintf(stderr, "\n");
1771 print_time(stderr, walltime_accounting, step, ir, cr);
1774 /* Ion/water position swapping.
1775 * Not done in last step since trajectory writing happens before this call
1776 * in the MD loop and exchanges would be lost anyway. */
1777 bNeedRepartition = FALSE;
1778 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1779 do_per_step(step, ir->swap->nstswap))
1781 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1782 bRerunMD ? rerun_fr.x : state->x,
1783 bRerunMD ? rerun_fr.box : state->box,
1784 top_global, MASTER(cr) && bVerbose, bRerunMD);
1786 if (bNeedRepartition && DOMAINDECOMP(cr))
1788 dd_collect_state(cr->dd, state, state_global);
1792 /* Replica exchange */
1796 bExchanged = replica_exchange(fplog, cr, repl_ex,
1797 state_global, enerd,
1801 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1803 dd_partition_system(fplog, step, cr, TRUE, 1,
1804 state_global, top_global, ir,
1805 state, &f, mdatoms, top, fr,
1806 vsite, shellfc, constr,
1807 nrnb, wcycle, FALSE);
1812 bStartingFromCpt = FALSE;
1814 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1815 /* With all integrators, except VV, we need to retain the pressure
1816 * at the current step for coupling at the next step.
1818 if ((state->flags & (1<<estPRES_PREV)) &&
1820 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1822 /* Store the pressure in t_state for pressure coupling
1823 * at the next MD step.
1825 copy_mat(pres, state->pres_prev);
1828 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1830 if ( (membed != NULL) && (!bLastStep) )
1832 rescale_membed(step_rel, membed, state_global->x);
1839 /* read next frame from input trajectory */
1840 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1845 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1849 if (!bRerunMD || !rerun_fr.bStep)
1851 /* increase the MD step number */
1856 cycles = wallcycle_stop(wcycle, ewcSTEP);
1857 if (DOMAINDECOMP(cr) && wcycle)
1859 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1862 if (bPMETuneRunning || bPMETuneTry)
1864 /* PME grid + cut-off optimization with GPUs or PME nodes */
1866 /* Count the total cycles over the last steps */
1867 cycles_pmes += cycles;
1869 /* We can only switch cut-off at NS steps */
1870 if (step % ir->nstlist == 0)
1872 /* PME grid + cut-off optimization with GPUs or PME nodes */
1875 if (DDMASTER(cr->dd))
1877 /* PME node load is too high, start tuning */
1878 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1880 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1882 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1884 bPMETuneTry = FALSE;
1887 if (bPMETuneRunning)
1889 /* init_step might not be a multiple of nstlist,
1890 * but the first cycle is always skipped anyhow.
1893 pme_load_balance(pme_loadbal, cr,
1894 (bVerbose && MASTER(cr)) ? stderr : NULL,
1896 ir, state, cycles_pmes,
1897 fr->ic, fr->nbv, &fr->pmedata,
1900 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1901 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1902 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1903 fr->rlist = fr->ic->rlist;
1904 fr->rlistlong = fr->ic->rlistlong;
1905 fr->rcoulomb = fr->ic->rcoulomb;
1906 fr->rvdw = fr->ic->rvdw;
1908 if (ir->eDispCorr != edispcNO)
1910 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1917 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1918 gs.set[eglsRESETCOUNTERS] != 0)
1920 /* Reset all the counters related to performance over the run */
1921 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1922 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1923 wcycle_set_reset_counters(wcycle, -1);
1924 if (!(cr->duty & DUTY_PME))
1926 /* Tell our PME node to reset its counters */
1927 gmx_pme_send_resetcounters(cr, step);
1929 /* Correct max_hours for the elapsed time */
1930 max_hours -= elapsed_time/(60.0*60.0);
1931 bResetCountersHalfMaxH = FALSE;
1932 gs.set[eglsRESETCOUNTERS] = 0;
1935 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1936 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1939 /* End of main MD loop */
1942 /* Stop measuring walltime */
1943 walltime_accounting_end(walltime_accounting);
1945 if (bRerunMD && MASTER(cr))
1950 if (!(cr->duty & DUTY_PME))
1952 /* Tell the PME only node to finish */
1953 gmx_pme_send_finish(cr);
1958 if (ir->nstcalcenergy > 0 && !bRerunMD)
1960 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1961 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1968 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1970 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)));
1971 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1974 if (pme_loadbal != NULL)
1976 pme_loadbal_done(pme_loadbal, cr, fplog,
1977 fr->nbv != NULL && fr->nbv->bUseGPU);
1980 if (shellfc && fplog)
1982 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1983 (nconverged*100.0)/step_rel);
1984 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1988 if (repl_ex_nst > 0 && MASTER(cr))
1990 print_replica_exchange_statistics(fplog, repl_ex);
1993 /* IMD cleanup, if bIMD is TRUE. */
1994 IMD_finalize(ir->bIMD, ir->imd);
1996 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);