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.
54 #include "md_support.h"
55 #include "md_logging.h"
68 #include "domdec_network.h"
69 #include "gromacs/gmxlib/topsort.h"
73 #include "compute_io.h"
74 #include "checkpoint.h"
75 #include "mtop_util.h"
76 #include "sighandler.h"
79 #include "pme_loadbal.h"
82 #include "types/nlistheuristics.h"
83 #include "types/iteratedconstraints.h"
84 #include "nbnxn_cuda_data_mgmt.h"
86 #include "gromacs/utility/gmxmpi.h"
87 #include "gromacs/fileio/confio.h"
88 #include "gromacs/fileio/trajectory_writing.h"
89 #include "gromacs/fileio/trnio.h"
90 #include "gromacs/fileio/trxio.h"
91 #include "gromacs/fileio/xtcio.h"
92 #include "gromacs/timing/wallcycle.h"
93 #include "gromacs/timing/walltime_accounting.h"
94 #include "gromacs/pulling/pull.h"
95 #include "gromacs/swap/swapcoords.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", walltime_accounting);
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,
149 gmx_walltime_accounting_t walltime_accounting)
151 gmx_mdoutf_t outf = NULL;
152 gmx_int64_t step, step_rel;
154 double t, t0, lam0[efptNR];
155 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
156 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
157 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
158 bBornRadii, bStartingFromCpt;
159 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
160 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
161 bForceUpdate = FALSE, bCPT;
162 gmx_bool bMasterState;
163 int force_flags, cglo_flags;
164 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
169 t_state *bufstate = NULL;
170 matrix *scale_tot, pcoupl_mu, M, ebox;
173 gmx_repl_ex_t repl_ex = NULL;
176 t_mdebin *mdebin = NULL;
177 t_state *state = NULL;
178 rvec *f_global = NULL;
179 gmx_enerdata_t *enerd;
181 gmx_global_stat_t gstat;
182 gmx_update_t upd = NULL;
183 t_graph *graph = NULL;
185 gmx_rng_t mcrng = NULL;
186 gmx_groups_t *groups;
187 gmx_ekindata_t *ekind, *ekind_save;
188 gmx_shellfc_t shellfc;
189 int count, nconverged = 0;
192 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;
211 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
212 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
213 gmx_iterate_t iterate;
214 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
215 simulation stops. If equal to zero, don't
216 communicate any more between multisims.*/
217 /* PME load balancing data for GPU kernels */
218 pme_load_balancing_t pme_loadbal = NULL;
220 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
223 /* Temporary addition for FAHCORE checkpointing */
227 /* Check for special mdrun options */
228 bRerunMD = (Flags & MD_RERUN);
229 bAppend = (Flags & MD_APPENDFILES);
230 if (Flags & MD_RESETCOUNTERSHALFWAY)
234 /* Signal to reset the counters half the simulation steps. */
235 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
237 /* Signal to reset the counters halfway the simulation time. */
238 bResetCountersHalfMaxH = (max_hours > 0);
241 /* md-vv uses averaged full step velocities for T-control
242 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
243 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
245 if (bVV) /* to store the initial velocities while computing virial */
247 snew(cbuf, top_global->natoms);
249 /* all the iteratative cases - only if there are constraints */
250 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
251 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
252 false in this step. The correct value, true or false,
253 is set at each step, as it depends on the frequency of temperature
254 and pressure control.*/
255 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
259 /* Since we don't know if the frames read are related in any way,
260 * rebuild the neighborlist at every step.
263 ir->nstcalcenergy = 1;
267 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
269 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
270 bGStatEveryStep = (nstglobalcomm == 1);
272 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
275 "To reduce the energy communication with nstlist = -1\n"
276 "the neighbor list validity should not be checked at every step,\n"
277 "this means that exact integration is not guaranteed.\n"
278 "The neighbor list validity is checked after:\n"
279 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
280 "In most cases this will result in exact integration.\n"
281 "This reduces the energy communication by a factor of 2 to 3.\n"
282 "If you want less energy communication, set nstlist > 3.\n\n");
287 ir->nstxout_compressed = 0;
289 groups = &top_global->groups;
292 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
293 &(state_global->fep_state), lam0,
294 nrnb, top_global, &upd,
295 nfile, fnm, &outf, &mdebin,
296 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
298 clear_mat(total_vir);
300 /* Energy terms and groups */
302 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
304 if (DOMAINDECOMP(cr))
310 snew(f, top_global->natoms);
313 /* Kinetic energy data */
315 init_ekindata(fplog, top_global, &(ir->opts), ekind);
316 /* needed for iteration of constraints */
318 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
319 /* Copy the cos acceleration to the groups struct */
320 ekind->cosacc.cos_accel = ir->cos_accel;
322 gstat = global_stat_init(ir);
325 /* Check for polarizable models and flexible constraints */
326 shellfc = init_shell_flexcon(fplog, fr->cutoff_scheme == ecutsVERLET,
327 top_global, n_flexible_constraints(constr),
328 (ir->bContinuation ||
329 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
330 NULL : state_global->x);
334 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
335 set_deform_reference_box(upd,
336 deform_init_init_step_tpx,
337 deform_init_box_tpx);
338 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
342 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
343 if ((io > 2000) && MASTER(cr))
346 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
351 if (DOMAINDECOMP(cr))
353 top = dd_init_local_top(top_global);
356 dd_init_local_state(cr->dd, state_global, state);
358 if (DDMASTER(cr->dd) && ir->nstfout)
360 snew(f_global, state_global->natoms);
365 top = gmx_mtop_generate_local_top(top_global, ir);
367 forcerec_set_excl_load(fr, top);
369 state = serial_init_local_state(state_global);
372 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
376 set_vsite_top(vsite, top, mdatoms, cr);
379 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
381 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
386 make_local_shells(cr, mdatoms, shellfc);
389 setup_bonded_threading(fr, &top->idef);
392 if (DOMAINDECOMP(cr))
394 /* Distribute the charge groups over the nodes from the master node */
395 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
396 state_global, top_global, ir,
397 state, &f, mdatoms, top, fr,
398 vsite, shellfc, constr,
399 nrnb, wcycle, FALSE);
403 update_mdatoms(mdatoms, state->lambda[efptMASS]);
405 if (opt2bSet("-cpi", nfile, fnm))
407 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
411 bStateFromCP = FALSE;
416 init_expanded_ensemble(bStateFromCP, ir, &mcrng, &state->dfhist);
423 /* Update mdebin with energy history if appending to output files */
424 if (Flags & MD_APPENDFILES)
426 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
430 /* We might have read an energy history from checkpoint,
431 * free the allocated memory and reset the counts.
433 done_energyhistory(&state_global->enerhist);
434 init_energyhistory(&state_global->enerhist);
437 /* Set the initial energy history in state by updating once */
438 update_energyhistory(&state_global->enerhist, mdebin);
441 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
443 /* Set the random state if we read a checkpoint file */
444 set_stochd_state(upd, state);
447 if (state->flags & (1<<estMC_RNG))
449 set_mc_state(mcrng, state);
452 /* Initialize constraints */
453 if (constr && !DOMAINDECOMP(cr))
455 set_constraints(constr, top, ir, mdatoms, cr);
460 /* We need to be sure replica exchange can only occur
461 * when the energies are current */
462 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
463 "repl_ex_nst", &repl_ex_nst);
464 /* This check needs to happen before inter-simulation
465 * signals are initialized, too */
467 if (repl_ex_nst > 0 && MASTER(cr))
469 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
470 repl_ex_nst, repl_ex_nex, repl_ex_seed);
473 /* PME tuning is only supported with GPUs or PME nodes and not with rerun or LJ-PME.
474 * With perturbed charges with soft-core we should not change the cut-off.
476 if ((Flags & MD_TUNEPME) &&
477 EEL_PME(fr->eeltype) &&
478 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
479 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
480 !bRerunMD && !EVDW_PME(fr->vdwtype))
482 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
484 if (cr->duty & DUTY_PME)
486 /* Start tuning right away, as we can't measure the load */
487 bPMETuneRunning = TRUE;
491 /* Separate PME nodes, we can measure the PP/PME load balance */
496 if (!ir->bContinuation && !bRerunMD)
498 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
500 /* Set the velocities of frozen particles to zero */
501 for (i = 0; i < mdatoms->homenr; i++)
503 for (m = 0; m < DIM; m++)
505 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
515 /* Constrain the initial coordinates and velocities */
516 do_constrain_first(fplog, constr, ir, mdatoms, state,
521 /* Construct the virtual sites for the initial configuration */
522 construct_vsites(vsite, state->x, ir->delta_t, NULL,
523 top->idef.iparams, top->idef.il,
524 fr->ePBC, fr->bMolPBC, cr, state->box);
530 /* set free energy calculation frequency as the minimum
531 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
532 nstfep = ir->fepvals->nstdhdl;
535 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
539 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
542 /* I'm assuming we need global communication the first time! MRS */
543 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
544 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
545 | (bVV ? CGLO_PRESSURE : 0)
546 | (bVV ? CGLO_CONSTRAINT : 0)
547 | (bRerunMD ? CGLO_RERUNMD : 0)
548 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
550 bSumEkinhOld = FALSE;
551 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
552 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
553 constr, NULL, FALSE, state->box,
554 top_global, &bSumEkinhOld, cglo_flags);
555 if (ir->eI == eiVVAK)
557 /* a second call to get the half step temperature initialized as well */
558 /* we do the same call as above, but turn the pressure off -- internally to
559 compute_globals, this is recognized as a velocity verlet half-step
560 kinetic energy calculation. This minimized excess variables, but
561 perhaps loses some logic?*/
563 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
564 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
565 constr, NULL, FALSE, state->box,
566 top_global, &bSumEkinhOld,
567 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
570 /* Calculate the initial half step temperature, and save the ekinh_old */
571 if (!(Flags & MD_STARTFROMCPT))
573 for (i = 0; (i < ir->opts.ngtc); i++)
575 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
580 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
581 and there is no previous step */
584 /* if using an iterative algorithm, we need to create a working directory for the state. */
587 bufstate = init_bufstate(state);
590 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
591 temperature control */
592 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
596 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
599 "RMS relative constraint deviation after constraining: %.2e\n",
600 constr_rmsd(constr, FALSE));
602 if (EI_STATE_VELOCITY(ir->eI))
604 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
608 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
609 " input trajectory '%s'\n\n",
610 *(top_global->name), opt2fn("-rerun", nfile, fnm));
613 fprintf(stderr, "Calculated time to finish depends on nsteps from "
614 "run input file,\nwhich may not correspond to the time "
615 "needed to process input trajectory.\n\n");
621 fprintf(stderr, "starting mdrun '%s'\n",
622 *(top_global->name));
625 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
629 sprintf(tbuf, "%s", "infinite");
631 if (ir->init_step > 0)
633 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
634 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
635 gmx_step_str(ir->init_step, sbuf2),
636 ir->init_step*ir->delta_t);
640 fprintf(stderr, "%s steps, %s ps.\n",
641 gmx_step_str(ir->nsteps, sbuf), tbuf);
644 fprintf(fplog, "\n");
647 walltime_accounting_start(walltime_accounting);
648 wallcycle_start(wcycle, ewcRUN);
649 print_start(fplog, cr, walltime_accounting, "mdrun");
651 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
653 chkpt_ret = fcCheckPointParallel( cr->nodeid,
657 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
662 /***********************************************************
666 ************************************************************/
668 /* if rerunMD then read coordinates and velocities from input trajectory */
671 if (getenv("GMX_FORCE_UPDATE"))
679 bNotLastFrame = read_first_frame(oenv, &status,
680 opt2fn("-rerun", nfile, fnm),
681 &rerun_fr, TRX_NEED_X | TRX_READ_V);
682 if (rerun_fr.natoms != top_global->natoms)
685 "Number of atoms in trajectory (%d) does not match the "
686 "run input file (%d)\n",
687 rerun_fr.natoms, top_global->natoms);
689 if (ir->ePBC != epbcNONE)
693 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);
695 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
697 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
704 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
707 if (ir->ePBC != epbcNONE)
709 /* Set the shift vectors.
710 * Necessary here when have a static box different from the tpr box.
712 calc_shifts(rerun_fr.box, fr->shift_vec);
716 /* loop over MD steps or if rerunMD to end of input trajectory */
718 /* Skip the first Nose-Hoover integration when we get the state from tpx */
719 bStateFromTPX = !bStateFromCP;
720 bInitStep = bFirstStep && (bStateFromTPX || bVV);
721 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 if (ir->nstlist == -1)
734 init_nlistheuristics(&nlh, bGStatEveryStep, step);
737 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
739 /* check how many steps are left in other sims */
740 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
744 /* and stop now if we should */
745 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
746 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
747 while (!bLastStep || (bRerunMD && bNotLastFrame))
750 wallcycle_start(wcycle, ewcSTEP);
756 step = rerun_fr.step;
757 step_rel = step - ir->init_step;
770 bLastStep = (step_rel == ir->nsteps);
771 t = t0 + step*ir->delta_t;
774 if (ir->efep != efepNO || ir->bSimTemp)
776 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
777 requiring different logic. */
779 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
780 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
781 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
782 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
783 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
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)
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 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1105 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1106 f, bUpdateDoLR, fr->f_twin, fcd,
1107 ekind, M, upd, bInitStep, etrtVELOCITY1,
1108 cr, nrnb, constr, &top->idef);
1110 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1112 gmx_iterate_init(&iterate, TRUE);
1114 /* for iterations, we save these vectors, as we will be self-consistently iterating
1117 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1119 /* save the state */
1120 if (iterate.bIterationActive)
1122 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1125 bFirstIterate = TRUE;
1126 while (bFirstIterate || iterate.bIterationActive)
1128 if (iterate.bIterationActive)
1130 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1131 if (bFirstIterate && bTrotter)
1133 /* The first time through, we need a decent first estimate
1134 of veta(t+dt) to compute the constraints. Do
1135 this by computing the box volume part of the
1136 trotter integration at this time. Nothing else
1137 should be changed by this routine here. If
1138 !(first time), we start with the previous value
1141 veta_save = state->veta;
1142 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1143 vetanew = state->veta;
1144 state->veta = veta_save;
1149 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1151 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1152 state, fr->bMolPBC, graph, f,
1153 &top->idef, shake_vir,
1154 cr, nrnb, wcycle, upd, constr,
1155 TRUE, bCalcVir, vetanew);
1159 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1165 /* Need to unshift here if a do_force has been
1166 called in the previous step */
1167 unshift_self(graph, state->box, state->x);
1170 /* if VV, compute the pressure and constraints */
1171 /* For VV2, we strictly only need this if using pressure
1172 * control, but we really would like to have accurate pressures
1174 * Think about ways around this in the future?
1175 * For now, keep this choice in comments.
1177 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1178 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1180 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1181 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1183 bSumEkinhOld = TRUE;
1185 /* for vv, the first half of the integration actually corresponds to the previous step.
1186 So we need information from the last step in the first half of the integration */
1187 if (bGStat || do_per_step(step-1, nstglobalcomm))
1189 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1190 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1191 constr, NULL, FALSE, state->box,
1192 top_global, &bSumEkinhOld,
1195 | (bTemp ? CGLO_TEMPERATURE : 0)
1196 | (bPres ? CGLO_PRESSURE : 0)
1197 | (bPres ? CGLO_CONSTRAINT : 0)
1198 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1199 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1202 /* explanation of above:
1203 a) We compute Ekin at the full time step
1204 if 1) we are using the AveVel Ekin, and it's not the
1205 initial step, or 2) if we are using AveEkin, but need the full
1206 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1207 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1208 EkinAveVel because it's needed for the pressure */
1210 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1215 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1216 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1223 /* We need the kinetic energy at minus the half step for determining
1224 * the full step kinetic energy and possibly for T-coupling.*/
1225 /* This may not be quite working correctly yet . . . . */
1226 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1227 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1228 constr, NULL, FALSE, state->box,
1229 top_global, &bSumEkinhOld,
1230 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1235 if (iterate.bIterationActive &&
1236 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1237 state->veta, &vetanew))
1241 bFirstIterate = FALSE;
1244 if (bTrotter && !bInitStep)
1246 copy_mat(shake_vir, state->svir_prev);
1247 copy_mat(force_vir, state->fvir_prev);
1248 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1250 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1251 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1252 enerd->term[F_EKIN] = trace(ekind->ekin);
1255 /* if it's the initial step, we performed this first step just to get the constraint virial */
1256 if (bInitStep && ir->eI == eiVV)
1258 copy_rvecn(cbuf, state->v, 0, state->natoms);
1262 /* MRS -- now done iterating -- compute the conserved quantity */
1265 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1268 last_ekin = enerd->term[F_EKIN];
1270 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1272 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1274 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1277 sum_dhdl(enerd, state->lambda, ir->fepvals);
1281 /* ######## END FIRST UPDATE STEP ############## */
1282 /* ######## If doing VV, we now have v(dt) ###### */
1285 /* perform extended ensemble sampling in lambda - we don't
1286 actually move to the new state before outputting
1287 statistics, but if performing simulated tempering, we
1288 do update the velocities and the tau_t. */
1290 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, mcrng, state->v, mdatoms);
1291 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1292 copy_df_history(&state_global->dfhist, &state->dfhist);
1295 /* Now we have the energies and forces corresponding to the
1296 * coordinates at time t. We must output all of this before
1299 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1300 ir, state, state_global, top_global, fr, upd,
1301 outf, mdebin, ekind, f, f_global,
1302 wcycle, mcrng, &nchkpt,
1303 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1306 /* kludge -- virial is lost with restart for NPT control. Must restart */
1307 if (bStartingFromCpt && bVV)
1309 copy_mat(state->svir_prev, shake_vir);
1310 copy_mat(state->fvir_prev, force_vir);
1313 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1315 /* Check whether everything is still allright */
1316 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1317 #ifdef GMX_THREAD_MPI
1322 /* this is just make gs.sig compatible with the hack
1323 of sending signals around by MPI_Reduce with together with
1325 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1327 gs.sig[eglsSTOPCOND] = 1;
1329 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1331 gs.sig[eglsSTOPCOND] = -1;
1333 /* < 0 means stop at next step, > 0 means stop at next NS step */
1337 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1338 gmx_get_signal_name(),
1339 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1343 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1344 gmx_get_signal_name(),
1345 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1347 handled_stop_condition = (int)gmx_get_stop_condition();
1349 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1350 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1351 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1353 /* Signal to terminate the run */
1354 gs.sig[eglsSTOPCOND] = 1;
1357 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1359 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1362 if (bResetCountersHalfMaxH && MASTER(cr) &&
1363 elapsed_time > max_hours*60.0*60.0*0.495)
1365 gs.sig[eglsRESETCOUNTERS] = 1;
1368 if (ir->nstlist == -1 && !bRerunMD)
1370 /* When bGStatEveryStep=FALSE, global_stat is only called
1371 * when we check the atom displacements, not at NS steps.
1372 * This means that also the bonded interaction count check is not
1373 * performed immediately after NS. Therefore a few MD steps could
1374 * be performed with missing interactions.
1375 * But wrong energies are never written to file,
1376 * since energies are only written after global_stat
1379 if (step >= nlh.step_nscheck)
1381 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1382 nlh.scale_tot, state->x);
1386 /* This is not necessarily true,
1387 * but step_nscheck is determined quite conservatively.
1393 /* In parallel we only have to check for checkpointing in steps
1394 * where we do global communication,
1395 * otherwise the other nodes don't know.
1397 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1400 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1401 gs.set[eglsCHKPT] == 0)
1403 gs.sig[eglsCHKPT] = 1;
1406 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1411 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1413 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1415 gmx_bool bIfRandomize;
1416 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1417 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1418 if (constr && bIfRandomize)
1420 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1421 state, fr->bMolPBC, graph, f,
1422 &top->idef, tmp_vir,
1423 cr, nrnb, wcycle, upd, constr,
1424 TRUE, bCalcVir, vetanew);
1429 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1431 gmx_iterate_init(&iterate, TRUE);
1432 /* for iterations, we save these vectors, as we will be redoing the calculations */
1433 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1436 bFirstIterate = TRUE;
1437 while (bFirstIterate || iterate.bIterationActive)
1439 /* We now restore these vectors to redo the calculation with improved extended variables */
1440 if (iterate.bIterationActive)
1442 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1445 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1446 so scroll down for that logic */
1448 /* ######### START SECOND UPDATE STEP ################# */
1449 /* Box is changed in update() when we do pressure coupling,
1450 * but we should still use the old box for energy corrections and when
1451 * writing it to the energy file, so it matches the trajectory files for
1452 * the same timestep above. Make a copy in a separate array.
1454 copy_mat(state->box, lastbox);
1459 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1461 wallcycle_start(wcycle, ewcUPDATE);
1462 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1465 if (iterate.bIterationActive)
1473 /* we use a new value of scalevir to converge the iterations faster */
1474 scalevir = tracevir/trace(shake_vir);
1476 msmul(shake_vir, scalevir, shake_vir);
1477 m_add(force_vir, shake_vir, total_vir);
1478 clear_mat(shake_vir);
1480 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1481 /* We can only do Berendsen coupling after we have summed
1482 * the kinetic energy or virial. Since the happens
1483 * in global_state after update, we should only do it at
1484 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1489 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1490 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1495 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1497 /* velocity half-step update */
1498 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1499 bUpdateDoLR, fr->f_twin, fcd,
1500 ekind, M, upd, FALSE, etrtVELOCITY2,
1501 cr, nrnb, constr, &top->idef);
1504 /* Above, initialize just copies ekinh into ekin,
1505 * it doesn't copy position (for VV),
1506 * and entire integrator for MD.
1509 if (ir->eI == eiVVAK)
1511 copy_rvecn(state->x, cbuf, 0, state->natoms);
1513 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1515 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1516 bUpdateDoLR, fr->f_twin, fcd,
1517 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1518 wallcycle_stop(wcycle, ewcUPDATE);
1520 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1521 fr->bMolPBC, graph, f,
1522 &top->idef, shake_vir,
1523 cr, nrnb, wcycle, upd, constr,
1524 FALSE, bCalcVir, state->veta);
1526 if (ir->eI == eiVVAK)
1528 /* erase F_EKIN and F_TEMP here? */
1529 /* just compute the kinetic energy at the half step to perform a trotter step */
1530 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1531 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1532 constr, NULL, FALSE, lastbox,
1533 top_global, &bSumEkinhOld,
1534 cglo_flags | CGLO_TEMPERATURE
1536 wallcycle_start(wcycle, ewcUPDATE);
1537 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1538 /* now we know the scaling, we can compute the positions again again */
1539 copy_rvecn(cbuf, state->x, 0, state->natoms);
1541 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1543 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1544 bUpdateDoLR, fr->f_twin, fcd,
1545 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1546 wallcycle_stop(wcycle, ewcUPDATE);
1548 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1549 /* are the small terms in the shake_vir here due
1550 * to numerical errors, or are they important
1551 * physically? I'm thinking they are just errors, but not completely sure.
1552 * For now, will call without actually constraining, constr=NULL*/
1553 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1554 state, fr->bMolPBC, graph, f,
1555 &top->idef, tmp_vir,
1556 cr, nrnb, wcycle, upd, NULL,
1562 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1565 if (fr->bSepDVDL && fplog && do_log)
1567 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1571 /* this factor or 2 correction is necessary
1572 because half of the constraint force is removed
1573 in the vv step, so we have to double it. See
1574 the Redmine issue #1255. It is not yet clear
1575 if the factor of 2 is exact, or just a very
1576 good approximation, and this will be
1577 investigated. The next step is to see if this
1578 can be done adding a dhdl contribution from the
1579 rattle step, but this is somewhat more
1580 complicated with the current code. Will be
1581 investigated, hopefully for 4.6.3. However,
1582 this current solution is much better than
1583 having it completely wrong.
1585 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1589 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1594 /* Need to unshift here */
1595 unshift_self(graph, state->box, state->x);
1600 wallcycle_start(wcycle, ewcVSITECONSTR);
1603 shift_self(graph, state->box, state->x);
1605 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1606 top->idef.iparams, top->idef.il,
1607 fr->ePBC, fr->bMolPBC, cr, state->box);
1611 unshift_self(graph, state->box, state->x);
1613 wallcycle_stop(wcycle, ewcVSITECONSTR);
1616 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1617 /* With Leap-Frog we can skip compute_globals at
1618 * non-communication steps, but we need to calculate
1619 * the kinetic energy one step before communication.
1621 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1623 if (ir->nstlist == -1 && bFirstIterate)
1625 gs.sig[eglsNABNSB] = nlh.nabnsb;
1627 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1628 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1630 bFirstIterate ? &gs : NULL,
1631 (step_rel % gs.nstms == 0) &&
1632 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1634 top_global, &bSumEkinhOld,
1636 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1637 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1638 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1639 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1640 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1641 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1644 if (ir->nstlist == -1 && bFirstIterate)
1646 nlh.nabnsb = gs.set[eglsNABNSB];
1647 gs.set[eglsNABNSB] = 0;
1650 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1651 /* ############# END CALC EKIN AND PRESSURE ################# */
1653 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1654 the virial that should probably be addressed eventually. state->veta has better properies,
1655 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1656 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1658 if (iterate.bIterationActive &&
1659 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1660 trace(shake_vir), &tracevir))
1664 bFirstIterate = FALSE;
1667 if (!bVV || bRerunMD)
1669 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1670 sum_dhdl(enerd, state->lambda, ir->fepvals);
1672 update_box(fplog, step, ir, mdatoms, state, f,
1673 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1675 /* ################# END UPDATE STEP 2 ################# */
1676 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1678 /* The coordinates (x) were unshifted in update */
1681 /* We will not sum ekinh_old,
1682 * so signal that we still have to do it.
1684 bSumEkinhOld = TRUE;
1687 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1689 /* use the directly determined last velocity, not actually the averaged half steps */
1690 if (bTrotter && ir->eI == eiVV)
1692 enerd->term[F_EKIN] = last_ekin;
1694 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1698 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1702 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1704 /* ######### END PREPARING EDR OUTPUT ########### */
1709 gmx_bool do_dr, do_or;
1711 if (fplog && do_log && bDoExpanded)
1713 /* only needed if doing expanded ensemble */
1714 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1715 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1717 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1721 upd_mdebin(mdebin, bDoDHDL, TRUE,
1722 t, mdatoms->tmass, enerd, state,
1723 ir->fepvals, ir->expandedvals, lastbox,
1724 shake_vir, force_vir, total_vir, pres,
1725 ekind, mu_tot, constr);
1729 upd_mdebin_step(mdebin);
1732 do_dr = do_per_step(step, ir->nstdisreout);
1733 do_or = do_per_step(step, ir->nstorireout);
1735 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1737 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1739 if (ir->ePull != epullNO)
1741 pull_print_output(ir->pull, step, t);
1744 if (do_per_step(step, ir->nstlog))
1746 if (fflush(fplog) != 0)
1748 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1754 /* Have to do this part _after_ outputting the logfile and the edr file */
1755 /* Gets written into the state at the beginning of next loop*/
1756 state->fep_state = lamnew;
1758 /* Print the remaining wall clock time for the run */
1759 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1763 fprintf(stderr, "\n");
1765 print_time(stderr, walltime_accounting, step, ir, cr);
1768 /* Ion/water position swapping.
1769 * Not done in last step since trajectory writing happens before this call
1770 * in the MD loop and exchanges would be lost anyway. */
1771 bNeedRepartition = FALSE;
1772 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1773 do_per_step(step, ir->swap->nstswap))
1775 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1776 bRerunMD ? rerun_fr.x : state->x,
1777 bRerunMD ? rerun_fr.box : state->box,
1778 top_global, MASTER(cr) && bVerbose, bRerunMD);
1780 if (bNeedRepartition && DOMAINDECOMP(cr))
1782 dd_collect_state(cr->dd, state, state_global);
1786 /* Replica exchange */
1788 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1789 do_per_step(step, repl_ex_nst))
1791 bExchanged = replica_exchange(fplog, cr, repl_ex,
1792 state_global, enerd,
1796 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1798 dd_partition_system(fplog, step, cr, TRUE, 1,
1799 state_global, top_global, ir,
1800 state, &f, mdatoms, top, fr,
1801 vsite, shellfc, constr,
1802 nrnb, wcycle, FALSE);
1807 bStartingFromCpt = FALSE;
1809 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1810 /* With all integrators, except VV, we need to retain the pressure
1811 * at the current step for coupling at the next step.
1813 if ((state->flags & (1<<estPRES_PREV)) &&
1815 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1817 /* Store the pressure in t_state for pressure coupling
1818 * at the next MD step.
1820 copy_mat(pres, state->pres_prev);
1823 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1825 if ( (membed != NULL) && (!bLastStep) )
1827 rescale_membed(step_rel, membed, state_global->x);
1834 /* read next frame from input trajectory */
1835 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1840 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1844 if (!bRerunMD || !rerun_fr.bStep)
1846 /* increase the MD step number */
1851 cycles = wallcycle_stop(wcycle, ewcSTEP);
1852 if (DOMAINDECOMP(cr) && wcycle)
1854 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1857 if (bPMETuneRunning || bPMETuneTry)
1859 /* PME grid + cut-off optimization with GPUs or PME nodes */
1861 /* Count the total cycles over the last steps */
1862 cycles_pmes += cycles;
1864 /* We can only switch cut-off at NS steps */
1865 if (step % ir->nstlist == 0)
1867 /* PME grid + cut-off optimization with GPUs or PME nodes */
1870 if (DDMASTER(cr->dd))
1872 /* PME node load is too high, start tuning */
1873 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1875 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1877 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1879 bPMETuneTry = FALSE;
1882 if (bPMETuneRunning)
1884 /* init_step might not be a multiple of nstlist,
1885 * but the first cycle is always skipped anyhow.
1888 pme_load_balance(pme_loadbal, cr,
1889 (bVerbose && MASTER(cr)) ? stderr : NULL,
1891 ir, state, cycles_pmes,
1892 fr->ic, fr->nbv, &fr->pmedata,
1895 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1896 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1897 fr->rlist = fr->ic->rlist;
1898 fr->rlistlong = fr->ic->rlistlong;
1899 fr->rcoulomb = fr->ic->rcoulomb;
1900 fr->rvdw = fr->ic->rvdw;
1906 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1907 gs.set[eglsRESETCOUNTERS] != 0)
1909 /* Reset all the counters related to performance over the run */
1910 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1911 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1912 wcycle_set_reset_counters(wcycle, -1);
1913 if (!(cr->duty & DUTY_PME))
1915 /* Tell our PME node to reset its counters */
1916 gmx_pme_send_resetcounters(cr, step);
1918 /* Correct max_hours for the elapsed time */
1919 max_hours -= elapsed_time/(60.0*60.0);
1920 bResetCountersHalfMaxH = FALSE;
1921 gs.set[eglsRESETCOUNTERS] = 0;
1925 /* End of main MD loop */
1928 /* Stop measuring walltime */
1929 walltime_accounting_end(walltime_accounting);
1931 if (bRerunMD && MASTER(cr))
1936 if (!(cr->duty & DUTY_PME))
1938 /* Tell the PME only node to finish */
1939 gmx_pme_send_finish(cr);
1944 if (ir->nstcalcenergy > 0 && !bRerunMD)
1946 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1947 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1954 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1956 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)));
1957 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1960 if (pme_loadbal != NULL)
1962 pme_loadbal_done(pme_loadbal, cr, fplog,
1963 fr->nbv != NULL && fr->nbv->bUseGPU);
1966 if (shellfc && fplog)
1968 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1969 (nconverged*100.0)/step_rel);
1970 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1974 if (repl_ex_nst > 0 && MASTER(cr))
1976 print_replica_exchange_statistics(fplog, repl_ex);
1979 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);