2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
44 #include "gromacs/math/vec.h"
53 #include "md_support.h"
54 #include "md_logging.h"
66 #include "domdec_network.h"
70 #include "gromacs/gmxpreprocess/compute_io.h"
71 #include "checkpoint.h"
72 #include "gromacs/topology/mtop_util.h"
73 #include "sighandler.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "pme_loadbal.h"
79 #include "types/nlistheuristics.h"
80 #include "types/iteratedconstraints.h"
81 #include "nbnxn_cuda_data_mgmt.h"
83 #include "gromacs/fileio/confio.h"
84 #include "gromacs/fileio/mdoutf.h"
85 #include "gromacs/fileio/trajectory_writing.h"
86 #include "gromacs/fileio/trnio.h"
87 #include "gromacs/fileio/trxio.h"
88 #include "gromacs/fileio/xtcio.h"
89 #include "gromacs/imd/imd.h"
90 #include "gromacs/pbcutil/mshift.h"
91 #include "gromacs/pbcutil/pbc.h"
92 #include "gromacs/pulling/pull.h"
93 #include "gromacs/swap/swapcoords.h"
94 #include "gromacs/timing/wallcycle.h"
95 #include "gromacs/timing/walltime_accounting.h"
96 #include "gromacs/utility/gmxmpi.h"
97 #include "gromacs/utility/smalloc.h"
100 #include "corewrap.h"
103 static void reset_all_counters(FILE *fplog, t_commrec *cr,
105 gmx_int64_t *step_rel, t_inputrec *ir,
106 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
107 gmx_walltime_accounting_t walltime_accounting,
108 nbnxn_cuda_ptr_t cu_nbv)
110 char sbuf[STEPSTRSIZE];
112 /* Reset all the counters related to performance over the run */
113 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
114 gmx_step_str(step, sbuf));
118 nbnxn_cuda_reset_timings(cu_nbv);
121 wallcycle_stop(wcycle, ewcRUN);
122 wallcycle_reset_all(wcycle);
123 if (DOMAINDECOMP(cr))
125 reset_dd_statistics_counters(cr->dd);
128 ir->init_step += *step_rel;
129 ir->nsteps -= *step_rel;
131 wallcycle_start(wcycle, ewcRUN);
132 walltime_accounting_start(walltime_accounting);
133 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
136 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
137 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
139 gmx_vsite_t *vsite, gmx_constr_t constr,
140 int stepout, t_inputrec *ir,
141 gmx_mtop_t *top_global,
143 t_state *state_global,
145 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
146 gmx_edsam_t ed, t_forcerec *fr,
147 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
148 real cpt_period, real max_hours,
149 const char gmx_unused *deviceOptions,
152 gmx_walltime_accounting_t walltime_accounting)
154 gmx_mdoutf_t outf = NULL;
155 gmx_int64_t step, step_rel;
157 double t, t0, lam0[efptNR];
158 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
159 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
160 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
161 bBornRadii, bStartingFromCpt;
162 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
163 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
164 bForceUpdate = FALSE, bCPT;
165 gmx_bool bMasterState;
166 int force_flags, cglo_flags;
167 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
172 t_state *bufstate = NULL;
176 gmx_repl_ex_t repl_ex = NULL;
179 t_mdebin *mdebin = NULL;
180 t_state *state = NULL;
181 rvec *f_global = NULL;
182 gmx_enerdata_t *enerd;
184 gmx_global_stat_t gstat;
185 gmx_update_t upd = NULL;
186 t_graph *graph = NULL;
188 gmx_groups_t *groups;
189 gmx_ekindata_t *ekind, *ekind_save;
190 gmx_shellfc_t shellfc;
191 int count, nconverged = 0;
193 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged, bNeedRepartition;
194 gmx_bool bResetCountersHalfMaxH = FALSE;
195 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
196 gmx_bool bUpdateDoLR;
200 real veta_save, scalevir, tracevir;
206 real saved_conserved_quantity = 0;
210 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
211 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
212 gmx_iterate_t iterate;
213 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
214 simulation stops. If equal to zero, don't
215 communicate any more between multisims.*/
216 /* PME load balancing data for GPU kernels */
217 pme_load_balancing_t pme_loadbal = NULL;
219 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
222 gmx_bool bIMDstep = FALSE;
225 /* Temporary addition for FAHCORE checkpointing */
229 /* Check for special mdrun options */
230 bRerunMD = (Flags & MD_RERUN);
231 if (Flags & MD_RESETCOUNTERSHALFWAY)
235 /* Signal to reset the counters half the simulation steps. */
236 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
238 /* Signal to reset the counters halfway the simulation time. */
239 bResetCountersHalfMaxH = (max_hours > 0);
242 /* md-vv uses averaged full step velocities for T-control
243 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
244 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
246 if (bVV) /* to store the initial velocities while computing virial */
248 snew(cbuf, top_global->natoms);
250 /* all the iteratative cases - only if there are constraints */
251 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
252 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
253 false in this step. The correct value, true or false,
254 is set at each step, as it depends on the frequency of temperature
255 and pressure control.*/
256 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
260 /* Since we don't know if the frames read are related in any way,
261 * rebuild the neighborlist at every step.
264 ir->nstcalcenergy = 1;
268 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
270 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
271 bGStatEveryStep = (nstglobalcomm == 1);
273 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
276 "To reduce the energy communication with nstlist = -1\n"
277 "the neighbor list validity should not be checked at every step,\n"
278 "this means that exact integration is not guaranteed.\n"
279 "The neighbor list validity is checked after:\n"
280 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
281 "In most cases this will result in exact integration.\n"
282 "This reduces the energy communication by a factor of 2 to 3.\n"
283 "If you want less energy communication, set nstlist > 3.\n\n");
288 ir->nstxout_compressed = 0;
290 groups = &top_global->groups;
293 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
294 &(state_global->fep_state), lam0,
295 nrnb, top_global, &upd,
296 nfile, fnm, &outf, &mdebin,
297 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
299 clear_mat(total_vir);
301 /* Energy terms and groups */
303 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
305 if (DOMAINDECOMP(cr))
311 snew(f, top_global->natoms);
314 /* Kinetic energy data */
316 init_ekindata(fplog, top_global, &(ir->opts), ekind);
317 /* needed for iteration of constraints */
319 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
320 /* Copy the cos acceleration to the groups struct */
321 ekind->cosacc.cos_accel = ir->cos_accel;
323 gstat = global_stat_init(ir);
326 /* Check for polarizable models and flexible constraints */
327 shellfc = init_shell_flexcon(fplog, fr->cutoff_scheme == ecutsVERLET,
328 top_global, n_flexible_constraints(constr),
329 (ir->bContinuation ||
330 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
331 NULL : state_global->x);
333 if (shellfc && ir->eI == eiNM)
335 /* Currently shells don't work with Normal Modes */
336 gmx_fatal(FARGS, "Normal Mode analysis is not supported with shells.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
339 if (vsite && ir->eI == eiNM)
341 /* Currently virtual sites don't work with Normal Modes */
342 gmx_fatal(FARGS, "Normal Mode analysis is not supported with virtual sites.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
347 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
348 set_deform_reference_box(upd,
349 deform_init_init_step_tpx,
350 deform_init_box_tpx);
351 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
355 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
356 if ((io > 2000) && MASTER(cr))
359 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
364 if (DOMAINDECOMP(cr))
366 top = dd_init_local_top(top_global);
369 dd_init_local_state(cr->dd, state_global, state);
371 if (DDMASTER(cr->dd) && ir->nstfout)
373 snew(f_global, state_global->natoms);
378 top = gmx_mtop_generate_local_top(top_global, ir);
380 forcerec_set_excl_load(fr, top);
382 state = serial_init_local_state(state_global);
385 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
389 set_vsite_top(vsite, top, mdatoms, cr);
392 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
394 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
399 make_local_shells(cr, mdatoms, shellfc);
402 setup_bonded_threading(fr, &top->idef);
405 /* Set up interactive MD (IMD) */
406 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
407 nfile, fnm, oenv, imdport, Flags);
409 if (DOMAINDECOMP(cr))
411 /* Distribute the charge groups over the nodes from the master node */
412 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
413 state_global, top_global, ir,
414 state, &f, mdatoms, top, fr,
415 vsite, shellfc, constr,
416 nrnb, wcycle, FALSE);
420 update_mdatoms(mdatoms, state->lambda[efptMASS]);
422 if (opt2bSet("-cpi", nfile, fnm))
424 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
428 bStateFromCP = FALSE;
433 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
440 /* Update mdebin with energy history if appending to output files */
441 if (Flags & MD_APPENDFILES)
443 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
447 /* We might have read an energy history from checkpoint,
448 * free the allocated memory and reset the counts.
450 done_energyhistory(&state_global->enerhist);
451 init_energyhistory(&state_global->enerhist);
454 /* Set the initial energy history in state by updating once */
455 update_energyhistory(&state_global->enerhist, mdebin);
458 /* Initialize constraints */
459 if (constr && !DOMAINDECOMP(cr))
461 set_constraints(constr, top, ir, mdatoms, cr);
466 /* We need to be sure replica exchange can only occur
467 * when the energies are current */
468 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
469 "repl_ex_nst", &repl_ex_nst);
470 /* This check needs to happen before inter-simulation
471 * signals are initialized, too */
473 if (repl_ex_nst > 0 && MASTER(cr))
475 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
476 repl_ex_nst, repl_ex_nex, repl_ex_seed);
479 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
480 * PME tuning is not supported with PME only for LJ and not for Coulomb.
482 if ((Flags & MD_TUNEPME) &&
483 EEL_PME(fr->eeltype) &&
484 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
487 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
489 if (cr->duty & DUTY_PME)
491 /* Start tuning right away, as we can't measure the load */
492 bPMETuneRunning = TRUE;
496 /* Separate PME nodes, we can measure the PP/PME load balance */
501 if (!ir->bContinuation && !bRerunMD)
503 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
505 /* Set the velocities of frozen particles to zero */
506 for (i = 0; i < mdatoms->homenr; i++)
508 for (m = 0; m < DIM; m++)
510 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
520 /* Constrain the initial coordinates and velocities */
521 do_constrain_first(fplog, constr, ir, mdatoms, state,
526 /* Construct the virtual sites for the initial configuration */
527 construct_vsites(vsite, state->x, ir->delta_t, NULL,
528 top->idef.iparams, top->idef.il,
529 fr->ePBC, fr->bMolPBC, cr, state->box);
535 /* set free energy calculation frequency as the minimum
536 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
537 nstfep = ir->fepvals->nstdhdl;
540 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
544 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
547 /* I'm assuming we need global communication the first time! MRS */
548 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
549 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
550 | (bVV ? CGLO_PRESSURE : 0)
551 | (bVV ? CGLO_CONSTRAINT : 0)
552 | (bRerunMD ? CGLO_RERUNMD : 0)
553 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
555 bSumEkinhOld = FALSE;
556 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
557 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
558 constr, NULL, FALSE, state->box,
559 top_global, &bSumEkinhOld, cglo_flags);
560 if (ir->eI == eiVVAK)
562 /* a second call to get the half step temperature initialized as well */
563 /* we do the same call as above, but turn the pressure off -- internally to
564 compute_globals, this is recognized as a velocity verlet half-step
565 kinetic energy calculation. This minimized excess variables, but
566 perhaps loses some logic?*/
568 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
569 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
570 constr, NULL, FALSE, state->box,
571 top_global, &bSumEkinhOld,
572 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
575 /* Calculate the initial half step temperature, and save the ekinh_old */
576 if (!(Flags & MD_STARTFROMCPT))
578 for (i = 0; (i < ir->opts.ngtc); i++)
580 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
585 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
586 and there is no previous step */
589 /* if using an iterative algorithm, we need to create a working directory for the state. */
592 bufstate = init_bufstate(state);
595 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
596 temperature control */
597 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
601 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
604 "RMS relative constraint deviation after constraining: %.2e\n",
605 constr_rmsd(constr, FALSE));
607 if (EI_STATE_VELOCITY(ir->eI))
609 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
613 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
614 " input trajectory '%s'\n\n",
615 *(top_global->name), opt2fn("-rerun", nfile, fnm));
618 fprintf(stderr, "Calculated time to finish depends on nsteps from "
619 "run input file,\nwhich may not correspond to the time "
620 "needed to process input trajectory.\n\n");
626 fprintf(stderr, "starting mdrun '%s'\n",
627 *(top_global->name));
630 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
634 sprintf(tbuf, "%s", "infinite");
636 if (ir->init_step > 0)
638 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
639 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
640 gmx_step_str(ir->init_step, sbuf2),
641 ir->init_step*ir->delta_t);
645 fprintf(stderr, "%s steps, %s ps.\n",
646 gmx_step_str(ir->nsteps, sbuf), tbuf);
649 fprintf(fplog, "\n");
652 walltime_accounting_start(walltime_accounting);
653 wallcycle_start(wcycle, ewcRUN);
654 print_start(fplog, cr, walltime_accounting, "mdrun");
656 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
658 chkpt_ret = fcCheckPointParallel( cr->nodeid,
662 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
667 /***********************************************************
671 ************************************************************/
673 /* if rerunMD then read coordinates and velocities from input trajectory */
676 if (getenv("GMX_FORCE_UPDATE"))
684 bNotLastFrame = read_first_frame(oenv, &status,
685 opt2fn("-rerun", nfile, fnm),
686 &rerun_fr, TRX_NEED_X | TRX_READ_V);
687 if (rerun_fr.natoms != top_global->natoms)
690 "Number of atoms in trajectory (%d) does not match the "
691 "run input file (%d)\n",
692 rerun_fr.natoms, top_global->natoms);
694 if (ir->ePBC != epbcNONE)
698 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);
700 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
702 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
709 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
712 if (ir->ePBC != epbcNONE)
714 /* Set the shift vectors.
715 * Necessary here when have a static box different from the tpr box.
717 calc_shifts(rerun_fr.box, fr->shift_vec);
721 /* loop over MD steps or if rerunMD to end of input trajectory */
723 /* Skip the first Nose-Hoover integration when we get the state from tpx */
724 bStateFromTPX = !bStateFromCP;
725 bInitStep = bFirstStep && (bStateFromTPX || bVV);
726 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
727 bSumEkinhOld = FALSE;
729 bNeedRepartition = FALSE;
731 init_global_signals(&gs, cr, ir, repl_ex_nst);
733 step = ir->init_step;
736 init_nlistheuristics(&nlh, bGStatEveryStep, step);
738 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
740 /* check how many steps are left in other sims */
741 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
745 /* and stop now if we should */
746 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
747 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
748 while (!bLastStep || (bRerunMD && bNotLastFrame))
751 wallcycle_start(wcycle, ewcSTEP);
757 step = rerun_fr.step;
758 step_rel = step - ir->init_step;
771 bLastStep = (step_rel == ir->nsteps);
772 t = t0 + step*ir->delta_t;
775 if (ir->efep != efepNO || ir->bSimTemp)
777 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
778 requiring different logic. */
780 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
781 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
782 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
783 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
784 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
789 update_annealing_target_temp(&(ir->opts), t);
794 if (!DOMAINDECOMP(cr) || MASTER(cr))
796 for (i = 0; i < state_global->natoms; i++)
798 copy_rvec(rerun_fr.x[i], state_global->x[i]);
802 for (i = 0; i < state_global->natoms; i++)
804 copy_rvec(rerun_fr.v[i], state_global->v[i]);
809 for (i = 0; i < state_global->natoms; i++)
811 clear_rvec(state_global->v[i]);
815 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
816 " Ekin, temperature and pressure are incorrect,\n"
817 " the virial will be incorrect when constraints are present.\n"
819 bRerunWarnNoV = FALSE;
823 copy_mat(rerun_fr.box, state_global->box);
824 copy_mat(state_global->box, state->box);
826 if (vsite && (Flags & MD_RERUN_VSITE))
828 if (DOMAINDECOMP(cr))
830 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
834 /* Following is necessary because the graph may get out of sync
835 * with the coordinates if we only have every N'th coordinate set
837 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
838 shift_self(graph, state->box, state->x);
840 construct_vsites(vsite, state->x, ir->delta_t, state->v,
841 top->idef.iparams, top->idef.il,
842 fr->ePBC, fr->bMolPBC, cr, state->box);
845 unshift_self(graph, state->box, state->x);
850 /* Stop Center of Mass motion */
851 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
855 /* for rerun MD always do Neighbour Searching */
856 bNS = (bFirstStep || ir->nstlist != 0);
861 /* Determine whether or not to do Neighbour Searching and LR */
862 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
864 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
865 (ir->nstlist == -1 && nlh.nabnsb > 0));
867 if (bNS && ir->nstlist == -1)
869 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
873 /* check whether we should stop because another simulation has
877 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
878 (multisim_nsteps != ir->nsteps) )
885 "Stopping simulation %d because another one has finished\n",
889 gs.sig[eglsCHKPT] = 1;
894 /* < 0 means stop at next step, > 0 means stop at next NS step */
895 if ( (gs.set[eglsSTOPCOND] < 0) ||
896 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
901 /* Determine whether or not to update the Born radii if doing GB */
902 bBornRadii = bFirstStep;
903 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
908 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
909 do_verbose = bVerbose &&
910 (step % stepout == 0 || bFirstStep || bLastStep);
912 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
920 bMasterState = FALSE;
921 /* Correct the new box if it is too skewed */
922 if (DYNAMIC_BOX(*ir))
924 if (correct_box(fplog, step, state->box, graph))
929 if (DOMAINDECOMP(cr) && bMasterState)
931 dd_collect_state(cr->dd, state, state_global);
935 if (DOMAINDECOMP(cr))
937 /* Repartition the domain decomposition */
938 wallcycle_start(wcycle, ewcDOMDEC);
939 dd_partition_system(fplog, step, cr,
940 bMasterState, nstglobalcomm,
941 state_global, top_global, ir,
942 state, &f, mdatoms, top, fr,
943 vsite, shellfc, constr,
945 do_verbose && !bPMETuneRunning);
946 wallcycle_stop(wcycle, ewcDOMDEC);
947 /* If using an iterative integrator, reallocate space to match the decomposition */
951 if (MASTER(cr) && do_log)
953 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
956 if (ir->efep != efepNO)
958 update_mdatoms(mdatoms, state->lambda[efptMASS]);
961 if ((bRerunMD && rerun_fr.bV) || bExchanged)
964 /* We need the kinetic energy at minus the half step for determining
965 * the full step kinetic energy and possibly for T-coupling.*/
966 /* This may not be quite working correctly yet . . . . */
967 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
968 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
969 constr, NULL, FALSE, state->box,
970 top_global, &bSumEkinhOld,
971 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
973 clear_mat(force_vir);
975 /* We write a checkpoint at this MD step when:
976 * either at an NS step when we signalled through gs,
977 * or at the last step (but not when we do not want confout),
978 * but never at the first step or with rerun.
980 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
981 (bLastStep && (Flags & MD_CONFOUT))) &&
982 step > ir->init_step && !bRerunMD);
985 gs.set[eglsCHKPT] = 0;
988 /* Determine the energy and pressure:
989 * at nstcalcenergy steps and at energy output steps (set below).
991 if (EI_VV(ir->eI) && (!bInitStep))
993 /* for vv, the first half of the integration actually corresponds
994 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
995 but the virial needs to be calculated on both the current step and the 'next' step. Future
996 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
998 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
999 bCalcVir = bCalcEner ||
1000 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1004 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1005 bCalcVir = bCalcEner ||
1006 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1009 /* Do we need global communication ? */
1010 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1011 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1012 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1014 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1016 if (do_ene || do_log)
1023 /* these CGLO_ options remain the same throughout the iteration */
1024 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1025 (bGStat ? CGLO_GSTAT : 0)
1028 force_flags = (GMX_FORCE_STATECHANGED |
1029 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1030 GMX_FORCE_ALLFORCES |
1032 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1033 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1034 (bDoFEP ? GMX_FORCE_DHDL : 0)
1039 if (do_per_step(step, ir->nstcalclr))
1041 force_flags |= GMX_FORCE_DO_LR;
1047 /* Now is the time to relax the shells */
1048 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1049 ir, bNS, force_flags,
1052 state, f, force_vir, mdatoms,
1053 nrnb, wcycle, graph, groups,
1054 shellfc, fr, bBornRadii, t, mu_tot,
1056 mdoutf_get_fp_field(outf));
1066 /* The coordinates (x) are shifted (to get whole molecules)
1068 * This is parallellized as well, and does communication too.
1069 * Check comments in sim_util.c
1071 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1072 state->box, state->x, &state->hist,
1073 f, force_vir, mdatoms, enerd, fcd,
1074 state->lambda, graph,
1075 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1076 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1079 if (bVV && !bStartingFromCpt && !bRerunMD)
1080 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1082 if (ir->eI == eiVV && bInitStep)
1084 /* if using velocity verlet with full time step Ekin,
1085 * take the first half step only to compute the
1086 * virial for the first step. From there,
1087 * revert back to the initial coordinates
1088 * so that the input is actually the initial step.
1090 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1094 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1095 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1098 /* If we are using twin-range interactions where the long-range component
1099 * is only evaluated every nstcalclr>1 steps, we should do a special update
1100 * step to combine the long-range forces on these steps.
1101 * For nstcalclr=1 this is not done, since the forces would have been added
1102 * directly to the short-range forces already.
1104 * TODO Remove various aspects of VV+twin-range in master
1105 * branch, because VV integrators did not ever support
1106 * twin-range multiple time stepping with constraints.
1108 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1110 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1111 f, bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1112 ekind, M, upd, bInitStep, etrtVELOCITY1,
1113 cr, nrnb, constr, &top->idef);
1115 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1117 gmx_iterate_init(&iterate, TRUE);
1119 /* for iterations, we save these vectors, as we will be self-consistently iterating
1122 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1124 /* save the state */
1125 if (iterate.bIterationActive)
1127 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1130 bFirstIterate = TRUE;
1131 while (bFirstIterate || iterate.bIterationActive)
1133 if (iterate.bIterationActive)
1135 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1136 if (bFirstIterate && bTrotter)
1138 /* The first time through, we need a decent first estimate
1139 of veta(t+dt) to compute the constraints. Do
1140 this by computing the box volume part of the
1141 trotter integration at this time. Nothing else
1142 should be changed by this routine here. If
1143 !(first time), we start with the previous value
1146 veta_save = state->veta;
1147 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1148 vetanew = state->veta;
1149 state->veta = veta_save;
1154 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1156 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1157 state, fr->bMolPBC, graph, f,
1158 &top->idef, shake_vir,
1159 cr, nrnb, wcycle, upd, constr,
1160 TRUE, bCalcVir, vetanew);
1162 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1164 /* Correct the virial for multiple time stepping */
1165 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1170 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1176 /* Need to unshift here if a do_force has been
1177 called in the previous step */
1178 unshift_self(graph, state->box, state->x);
1181 /* if VV, compute the pressure and constraints */
1182 /* For VV2, we strictly only need this if using pressure
1183 * control, but we really would like to have accurate pressures
1185 * Think about ways around this in the future?
1186 * For now, keep this choice in comments.
1188 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1189 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1191 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1192 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1194 bSumEkinhOld = TRUE;
1196 /* for vv, the first half of the integration actually corresponds to the previous step.
1197 So we need information from the last step in the first half of the integration */
1198 if (bGStat || do_per_step(step-1, nstglobalcomm))
1200 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1201 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1202 constr, NULL, FALSE, state->box,
1203 top_global, &bSumEkinhOld,
1206 | (bTemp ? CGLO_TEMPERATURE : 0)
1207 | (bPres ? CGLO_PRESSURE : 0)
1208 | (bPres ? CGLO_CONSTRAINT : 0)
1209 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1210 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1213 /* explanation of above:
1214 a) We compute Ekin at the full time step
1215 if 1) we are using the AveVel Ekin, and it's not the
1216 initial step, or 2) if we are using AveEkin, but need the full
1217 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1218 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1219 EkinAveVel because it's needed for the pressure */
1221 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1226 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1227 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1234 /* We need the kinetic energy at minus the half step for determining
1235 * the full step kinetic energy and possibly for T-coupling.*/
1236 /* This may not be quite working correctly yet . . . . */
1237 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1238 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1239 constr, NULL, FALSE, state->box,
1240 top_global, &bSumEkinhOld,
1241 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1246 if (iterate.bIterationActive &&
1247 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1248 state->veta, &vetanew))
1252 bFirstIterate = FALSE;
1255 if (bTrotter && !bInitStep)
1257 copy_mat(shake_vir, state->svir_prev);
1258 copy_mat(force_vir, state->fvir_prev);
1259 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1261 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1262 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1263 enerd->term[F_EKIN] = trace(ekind->ekin);
1266 /* if it's the initial step, we performed this first step just to get the constraint virial */
1267 if (bInitStep && ir->eI == eiVV)
1269 copy_rvecn(cbuf, state->v, 0, state->natoms);
1273 /* MRS -- now done iterating -- compute the conserved quantity */
1276 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1279 last_ekin = enerd->term[F_EKIN];
1281 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1283 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1285 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1288 sum_dhdl(enerd, state->lambda, ir->fepvals);
1292 /* ######## END FIRST UPDATE STEP ############## */
1293 /* ######## If doing VV, we now have v(dt) ###### */
1296 /* perform extended ensemble sampling in lambda - we don't
1297 actually move to the new state before outputting
1298 statistics, but if performing simulated tempering, we
1299 do update the velocities and the tau_t. */
1301 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1302 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1303 copy_df_history(&state_global->dfhist, &state->dfhist);
1306 /* Now we have the energies and forces corresponding to the
1307 * coordinates at time t. We must output all of this before
1310 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1311 ir, state, state_global, top_global, fr,
1312 outf, mdebin, ekind, f, f_global,
1314 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1316 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1317 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1319 /* kludge -- virial is lost with restart for NPT control. Must restart */
1320 if (bStartingFromCpt && bVV)
1322 copy_mat(state->svir_prev, shake_vir);
1323 copy_mat(state->fvir_prev, force_vir);
1326 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1328 /* Check whether everything is still allright */
1329 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1330 #ifdef GMX_THREAD_MPI
1335 /* this is just make gs.sig compatible with the hack
1336 of sending signals around by MPI_Reduce with together with
1338 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1340 gs.sig[eglsSTOPCOND] = 1;
1342 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1344 gs.sig[eglsSTOPCOND] = -1;
1346 /* < 0 means stop at next step, > 0 means stop at next NS step */
1350 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1351 gmx_get_signal_name(),
1352 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1356 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1357 gmx_get_signal_name(),
1358 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1360 handled_stop_condition = (int)gmx_get_stop_condition();
1362 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1363 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1364 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1366 /* Signal to terminate the run */
1367 gs.sig[eglsSTOPCOND] = 1;
1370 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1372 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1375 if (bResetCountersHalfMaxH && MASTER(cr) &&
1376 elapsed_time > max_hours*60.0*60.0*0.495)
1378 gs.sig[eglsRESETCOUNTERS] = 1;
1381 if (ir->nstlist == -1 && !bRerunMD)
1383 /* When bGStatEveryStep=FALSE, global_stat is only called
1384 * when we check the atom displacements, not at NS steps.
1385 * This means that also the bonded interaction count check is not
1386 * performed immediately after NS. Therefore a few MD steps could
1387 * be performed with missing interactions.
1388 * But wrong energies are never written to file,
1389 * since energies are only written after global_stat
1392 if (step >= nlh.step_nscheck)
1394 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1395 nlh.scale_tot, state->x);
1399 /* This is not necessarily true,
1400 * but step_nscheck is determined quite conservatively.
1406 /* In parallel we only have to check for checkpointing in steps
1407 * where we do global communication,
1408 * otherwise the other nodes don't know.
1410 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1413 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1414 gs.set[eglsCHKPT] == 0)
1416 gs.sig[eglsCHKPT] = 1;
1419 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1424 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1426 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1428 gmx_bool bIfRandomize;
1429 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1430 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1431 if (constr && bIfRandomize)
1433 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1434 state, fr->bMolPBC, graph, f,
1435 &top->idef, tmp_vir,
1436 cr, nrnb, wcycle, upd, constr,
1437 TRUE, bCalcVir, vetanew);
1442 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1444 gmx_iterate_init(&iterate, TRUE);
1445 /* for iterations, we save these vectors, as we will be redoing the calculations */
1446 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1449 bFirstIterate = TRUE;
1450 while (bFirstIterate || iterate.bIterationActive)
1452 /* We now restore these vectors to redo the calculation with improved extended variables */
1453 if (iterate.bIterationActive)
1455 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1458 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1459 so scroll down for that logic */
1461 /* ######### START SECOND UPDATE STEP ################# */
1462 /* Box is changed in update() when we do pressure coupling,
1463 * but we should still use the old box for energy corrections and when
1464 * writing it to the energy file, so it matches the trajectory files for
1465 * the same timestep above. Make a copy in a separate array.
1467 copy_mat(state->box, lastbox);
1472 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1474 wallcycle_start(wcycle, ewcUPDATE);
1475 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1478 if (iterate.bIterationActive)
1486 /* we use a new value of scalevir to converge the iterations faster */
1487 scalevir = tracevir/trace(shake_vir);
1489 msmul(shake_vir, scalevir, shake_vir);
1490 m_add(force_vir, shake_vir, total_vir);
1491 clear_mat(shake_vir);
1493 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1494 /* We can only do Berendsen coupling after we have summed
1495 * the kinetic energy or virial. Since the happens
1496 * in global_state after update, we should only do it at
1497 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1502 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1503 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1508 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1510 /* velocity half-step update */
1511 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1512 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1513 ekind, M, upd, FALSE, etrtVELOCITY2,
1514 cr, nrnb, constr, &top->idef);
1517 /* Above, initialize just copies ekinh into ekin,
1518 * it doesn't copy position (for VV),
1519 * and entire integrator for MD.
1522 if (ir->eI == eiVVAK)
1524 copy_rvecn(state->x, cbuf, 0, state->natoms);
1526 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1528 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1529 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1530 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1531 wallcycle_stop(wcycle, ewcUPDATE);
1533 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1534 fr->bMolPBC, graph, f,
1535 &top->idef, shake_vir,
1536 cr, nrnb, wcycle, upd, constr,
1537 FALSE, bCalcVir, state->veta);
1539 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1541 /* Correct the virial for multiple time stepping */
1542 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1545 if (ir->eI == eiVVAK)
1547 /* erase F_EKIN and F_TEMP here? */
1548 /* just compute the kinetic energy at the half step to perform a trotter step */
1549 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1550 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1551 constr, NULL, FALSE, lastbox,
1552 top_global, &bSumEkinhOld,
1553 cglo_flags | CGLO_TEMPERATURE
1555 wallcycle_start(wcycle, ewcUPDATE);
1556 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1557 /* now we know the scaling, we can compute the positions again again */
1558 copy_rvecn(cbuf, state->x, 0, state->natoms);
1560 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1562 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1563 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1564 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1565 wallcycle_stop(wcycle, ewcUPDATE);
1567 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1568 /* are the small terms in the shake_vir here due
1569 * to numerical errors, or are they important
1570 * physically? I'm thinking they are just errors, but not completely sure.
1571 * For now, will call without actually constraining, constr=NULL*/
1572 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1573 state, fr->bMolPBC, graph, f,
1574 &top->idef, tmp_vir,
1575 cr, nrnb, wcycle, upd, NULL,
1581 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1584 if (fr->bSepDVDL && fplog && do_log)
1586 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1590 /* this factor or 2 correction is necessary
1591 because half of the constraint force is removed
1592 in the vv step, so we have to double it. See
1593 the Redmine issue #1255. It is not yet clear
1594 if the factor of 2 is exact, or just a very
1595 good approximation, and this will be
1596 investigated. The next step is to see if this
1597 can be done adding a dhdl contribution from the
1598 rattle step, but this is somewhat more
1599 complicated with the current code. Will be
1600 investigated, hopefully for 4.6.3. However,
1601 this current solution is much better than
1602 having it completely wrong.
1604 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1608 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1613 /* Need to unshift here */
1614 unshift_self(graph, state->box, state->x);
1619 wallcycle_start(wcycle, ewcVSITECONSTR);
1622 shift_self(graph, state->box, state->x);
1624 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1625 top->idef.iparams, top->idef.il,
1626 fr->ePBC, fr->bMolPBC, cr, state->box);
1630 unshift_self(graph, state->box, state->x);
1632 wallcycle_stop(wcycle, ewcVSITECONSTR);
1635 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1636 /* With Leap-Frog we can skip compute_globals at
1637 * non-communication steps, but we need to calculate
1638 * the kinetic energy one step before communication.
1640 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1642 if (ir->nstlist == -1 && bFirstIterate)
1644 gs.sig[eglsNABNSB] = nlh.nabnsb;
1646 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1647 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1649 bFirstIterate ? &gs : NULL,
1650 (step_rel % gs.nstms == 0) &&
1651 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1653 top_global, &bSumEkinhOld,
1655 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1656 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1657 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1658 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1659 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1660 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1663 if (ir->nstlist == -1 && bFirstIterate)
1665 nlh.nabnsb = gs.set[eglsNABNSB];
1666 gs.set[eglsNABNSB] = 0;
1669 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1670 /* ############# END CALC EKIN AND PRESSURE ################# */
1672 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1673 the virial that should probably be addressed eventually. state->veta has better properies,
1674 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1675 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1677 if (iterate.bIterationActive &&
1678 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1679 trace(shake_vir), &tracevir))
1683 bFirstIterate = FALSE;
1686 if (!bVV || bRerunMD)
1688 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1689 sum_dhdl(enerd, state->lambda, ir->fepvals);
1691 update_box(fplog, step, ir, mdatoms, state, f,
1692 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1694 /* ################# END UPDATE STEP 2 ################# */
1695 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1697 /* The coordinates (x) were unshifted in update */
1700 /* We will not sum ekinh_old,
1701 * so signal that we still have to do it.
1703 bSumEkinhOld = TRUE;
1706 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1708 /* use the directly determined last velocity, not actually the averaged half steps */
1709 if (bTrotter && ir->eI == eiVV)
1711 enerd->term[F_EKIN] = last_ekin;
1713 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1717 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1721 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1723 /* ######### END PREPARING EDR OUTPUT ########### */
1728 gmx_bool do_dr, do_or;
1730 if (fplog && do_log && bDoExpanded)
1732 /* only needed if doing expanded ensemble */
1733 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1734 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1736 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1740 upd_mdebin(mdebin, bDoDHDL, TRUE,
1741 t, mdatoms->tmass, enerd, state,
1742 ir->fepvals, ir->expandedvals, lastbox,
1743 shake_vir, force_vir, total_vir, pres,
1744 ekind, mu_tot, constr);
1748 upd_mdebin_step(mdebin);
1751 do_dr = do_per_step(step, ir->nstdisreout);
1752 do_or = do_per_step(step, ir->nstorireout);
1754 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1756 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1758 if (ir->ePull != epullNO)
1760 pull_print_output(ir->pull, step, t);
1763 if (do_per_step(step, ir->nstlog))
1765 if (fflush(fplog) != 0)
1767 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1773 /* Have to do this part _after_ outputting the logfile and the edr file */
1774 /* Gets written into the state at the beginning of next loop*/
1775 state->fep_state = lamnew;
1777 /* Print the remaining wall clock time for the run */
1778 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1782 fprintf(stderr, "\n");
1784 print_time(stderr, walltime_accounting, step, ir, cr);
1787 /* Ion/water position swapping.
1788 * Not done in last step since trajectory writing happens before this call
1789 * in the MD loop and exchanges would be lost anyway. */
1790 bNeedRepartition = FALSE;
1791 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1792 do_per_step(step, ir->swap->nstswap))
1794 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1795 bRerunMD ? rerun_fr.x : state->x,
1796 bRerunMD ? rerun_fr.box : state->box,
1797 top_global, MASTER(cr) && bVerbose, bRerunMD);
1799 if (bNeedRepartition && DOMAINDECOMP(cr))
1801 dd_collect_state(cr->dd, state, state_global);
1805 /* Replica exchange */
1807 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1808 do_per_step(step, repl_ex_nst))
1810 bExchanged = replica_exchange(fplog, cr, repl_ex,
1811 state_global, enerd,
1815 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1817 dd_partition_system(fplog, step, cr, TRUE, 1,
1818 state_global, top_global, ir,
1819 state, &f, mdatoms, top, fr,
1820 vsite, shellfc, constr,
1821 nrnb, wcycle, FALSE);
1826 bStartingFromCpt = FALSE;
1828 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1829 /* With all integrators, except VV, we need to retain the pressure
1830 * at the current step for coupling at the next step.
1832 if ((state->flags & (1<<estPRES_PREV)) &&
1834 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1836 /* Store the pressure in t_state for pressure coupling
1837 * at the next MD step.
1839 copy_mat(pres, state->pres_prev);
1842 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1844 if ( (membed != NULL) && (!bLastStep) )
1846 rescale_membed(step_rel, membed, state_global->x);
1853 /* read next frame from input trajectory */
1854 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1859 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1863 if (!bRerunMD || !rerun_fr.bStep)
1865 /* increase the MD step number */
1870 cycles = wallcycle_stop(wcycle, ewcSTEP);
1871 if (DOMAINDECOMP(cr) && wcycle)
1873 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1876 if (bPMETuneRunning || bPMETuneTry)
1878 /* PME grid + cut-off optimization with GPUs or PME nodes */
1880 /* Count the total cycles over the last steps */
1881 cycles_pmes += cycles;
1883 /* We can only switch cut-off at NS steps */
1884 if (step % ir->nstlist == 0)
1886 /* PME grid + cut-off optimization with GPUs or PME nodes */
1889 if (DDMASTER(cr->dd))
1891 /* PME node load is too high, start tuning */
1892 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1894 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1896 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1898 bPMETuneTry = FALSE;
1901 if (bPMETuneRunning)
1903 /* init_step might not be a multiple of nstlist,
1904 * but the first cycle is always skipped anyhow.
1907 pme_load_balance(pme_loadbal, cr,
1908 (bVerbose && MASTER(cr)) ? stderr : NULL,
1910 ir, state, cycles_pmes,
1911 fr->ic, fr->nbv, &fr->pmedata,
1914 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1915 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1916 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1917 fr->rlist = fr->ic->rlist;
1918 fr->rlistlong = fr->ic->rlistlong;
1919 fr->rcoulomb = fr->ic->rcoulomb;
1920 fr->rvdw = fr->ic->rvdw;
1922 if (ir->eDispCorr != edispcNO)
1924 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1931 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1932 gs.set[eglsRESETCOUNTERS] != 0)
1934 /* Reset all the counters related to performance over the run */
1935 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1936 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1937 wcycle_set_reset_counters(wcycle, -1);
1938 if (!(cr->duty & DUTY_PME))
1940 /* Tell our PME node to reset its counters */
1941 gmx_pme_send_resetcounters(cr, step);
1943 /* Correct max_hours for the elapsed time */
1944 max_hours -= elapsed_time/(60.0*60.0);
1945 bResetCountersHalfMaxH = FALSE;
1946 gs.set[eglsRESETCOUNTERS] = 0;
1949 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1950 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1953 /* End of main MD loop */
1956 /* Stop measuring walltime */
1957 walltime_accounting_end(walltime_accounting);
1959 if (bRerunMD && MASTER(cr))
1964 if (!(cr->duty & DUTY_PME))
1966 /* Tell the PME only node to finish */
1967 gmx_pme_send_finish(cr);
1972 if (ir->nstcalcenergy > 0 && !bRerunMD)
1974 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1975 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1982 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1984 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)));
1985 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1988 if (pme_loadbal != NULL)
1990 pme_loadbal_done(pme_loadbal, cr, fplog,
1991 fr->nbv != NULL && fr->nbv->bUseGPU);
1994 if (shellfc && fplog)
1996 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1997 (nconverged*100.0)/step_rel);
1998 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
2002 if (repl_ex_nst > 0 && MASTER(cr))
2004 print_replica_exchange_statistics(fplog, repl_ex);
2007 /* IMD cleanup, if bIMD is TRUE. */
2008 IMD_finalize(ir->bIMD, ir->imd);
2010 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);