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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
58 #include "md_support.h"
59 #include "md_logging.h"
74 #include "mpelogging.h"
76 #include "domdec_network.h"
82 #include "compute_io.h"
84 #include "checkpoint.h"
85 #include "mtop_util.h"
86 #include "sighandler.h"
89 #include "pme_loadbal.h"
92 #include "types/nlistheuristics.h"
93 #include "types/iteratedconstraints.h"
94 #include "nbnxn_cuda_data_mgmt.h"
104 #include "corewrap.h"
107 static void reset_all_counters(FILE *fplog, t_commrec *cr,
108 gmx_large_int_t step,
109 gmx_large_int_t *step_rel, t_inputrec *ir,
110 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
111 gmx_runtime_t *runtime,
112 nbnxn_cuda_ptr_t cu_nbv)
114 char sbuf[STEPSTRSIZE];
116 /* Reset all the counters related to performance over the run */
117 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
118 gmx_step_str(step, sbuf));
122 nbnxn_cuda_reset_timings(cu_nbv);
125 wallcycle_stop(wcycle, ewcRUN);
126 wallcycle_reset_all(wcycle);
127 if (DOMAINDECOMP(cr))
129 reset_dd_statistics_counters(cr->dd);
132 ir->init_step += *step_rel;
133 ir->nsteps -= *step_rel;
135 wallcycle_start(wcycle, ewcRUN);
136 runtime_start(runtime);
137 print_date_and_time(fplog, cr->nodeid, "Restarted time", runtime);
140 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
141 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
143 gmx_vsite_t *vsite, gmx_constr_t constr,
144 int stepout, t_inputrec *ir,
145 gmx_mtop_t *top_global,
147 t_state *state_global,
149 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
150 gmx_edsam_t ed, t_forcerec *fr,
151 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
152 real cpt_period, real max_hours,
153 const char *deviceOptions,
155 gmx_runtime_t *runtime)
158 gmx_large_int_t step, step_rel;
160 double t, t0, lam0[efptNR];
161 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
162 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
163 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
164 bBornRadii, bStartingFromCpt;
165 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
166 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
167 bForceUpdate = FALSE, bCPT;
169 gmx_bool bMasterState;
170 int force_flags, cglo_flags;
171 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
176 t_state *bufstate = NULL;
177 matrix *scale_tot, pcoupl_mu, M, ebox;
180 gmx_repl_ex_t repl_ex = NULL;
183 t_mdebin *mdebin = NULL;
184 t_state *state = NULL;
185 rvec *f_global = NULL;
188 gmx_enerdata_t *enerd;
190 gmx_global_stat_t gstat;
191 gmx_update_t upd = NULL;
192 t_graph *graph = NULL;
194 gmx_rng_t mcrng = NULL;
196 gmx_groups_t *groups;
197 gmx_ekindata_t *ekind, *ekind_save;
198 gmx_shellfc_t shellfc;
199 int count, nconverged = 0;
202 gmx_bool bIonize = FALSE;
203 gmx_bool bTCR = FALSE, bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
205 gmx_bool bResetCountersHalfMaxH = FALSE;
206 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
207 gmx_bool bUpdateDoLR;
208 real mu_aver = 0, dvdl_constr;
209 int a0, a1, gnx = 0, ii;
210 atom_id *grpindex = NULL;
212 t_coupl_rec *tcr = NULL;
213 rvec *xcopy = NULL, *vcopy = NULL, *cbuf = NULL;
214 matrix boxcopy = {{0}}, lastbox;
216 real fom, oldfom, veta_save, pcurr, scalevir, tracevir;
223 real saved_conserved_quantity = 0;
228 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
229 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
230 gmx_iterate_t iterate;
231 gmx_large_int_t multisim_nsteps = -1; /* number of steps to do before first multisim
232 simulation stops. If equal to zero, don't
233 communicate any more between multisims.*/
234 /* PME load balancing data for GPU kernels */
235 pme_load_balancing_t pme_loadbal = NULL;
237 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
240 /* Temporary addition for FAHCORE checkpointing */
244 /* Check for special mdrun options */
245 bRerunMD = (Flags & MD_RERUN);
246 bIonize = (Flags & MD_IONIZE);
247 bFFscan = (Flags & MD_FFSCAN);
248 bAppend = (Flags & MD_APPENDFILES);
249 if (Flags & MD_RESETCOUNTERSHALFWAY)
253 /* Signal to reset the counters half the simulation steps. */
254 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
256 /* Signal to reset the counters halfway the simulation time. */
257 bResetCountersHalfMaxH = (max_hours > 0);
260 /* md-vv uses averaged full step velocities for T-control
261 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
262 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
264 if (bVV) /* to store the initial velocities while computing virial */
266 snew(cbuf, top_global->natoms);
268 /* all the iteratative cases - only if there are constraints */
269 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
270 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
271 false in this step. The correct value, true or false,
272 is set at each step, as it depends on the frequency of temperature
273 and pressure control.*/
274 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
278 /* Since we don't know if the frames read are related in any way,
279 * rebuild the neighborlist at every step.
282 ir->nstcalcenergy = 1;
286 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
288 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
289 bGStatEveryStep = (nstglobalcomm == 1);
291 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
294 "To reduce the energy communication with nstlist = -1\n"
295 "the neighbor list validity should not be checked at every step,\n"
296 "this means that exact integration is not guaranteed.\n"
297 "The neighbor list validity is checked after:\n"
298 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
299 "In most cases this will result in exact integration.\n"
300 "This reduces the energy communication by a factor of 2 to 3.\n"
301 "If you want less energy communication, set nstlist > 3.\n\n");
304 if (bRerunMD || bFFscan)
308 groups = &top_global->groups;
311 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
312 &(state_global->fep_state), lam0,
313 nrnb, top_global, &upd,
314 nfile, fnm, &outf, &mdebin,
315 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, state_global, Flags);
317 clear_mat(total_vir);
319 /* Energy terms and groups */
321 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
323 if (DOMAINDECOMP(cr))
329 snew(f, top_global->natoms);
332 /* Kinetic energy data */
334 init_ekindata(fplog, top_global, &(ir->opts), ekind);
335 /* needed for iteration of constraints */
337 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
338 /* Copy the cos acceleration to the groups struct */
339 ekind->cosacc.cos_accel = ir->cos_accel;
341 gstat = global_stat_init(ir);
344 /* Check for polarizable models and flexible constraints */
345 shellfc = init_shell_flexcon(fplog, fr->cutoff_scheme == ecutsVERLET,
346 top_global, n_flexible_constraints(constr),
347 (ir->bContinuation ||
348 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
349 NULL : state_global->x);
353 #ifdef GMX_THREAD_MPI
354 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
356 set_deform_reference_box(upd,
357 deform_init_init_step_tpx,
358 deform_init_box_tpx);
359 #ifdef GMX_THREAD_MPI
360 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
365 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
366 if ((io > 2000) && MASTER(cr))
369 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
374 if (DOMAINDECOMP(cr))
376 top = dd_init_local_top(top_global);
379 dd_init_local_state(cr->dd, state_global, state);
381 if (DDMASTER(cr->dd) && ir->nstfout)
383 snew(f_global, state_global->natoms);
390 /* Initialize the particle decomposition and split the topology */
391 top = split_system(fplog, top_global, ir, cr);
393 pd_cg_range(cr, &fr->cg0, &fr->hcg);
394 pd_at_range(cr, &a0, &a1);
398 top = gmx_mtop_generate_local_top(top_global, ir);
401 a1 = top_global->natoms;
404 forcerec_set_excl_load(fr, top, cr);
406 state = partdec_init_local_state(cr, state_global);
409 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
413 set_vsite_top(vsite, top, mdatoms, cr);
416 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
418 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
423 make_local_shells(cr, mdatoms, shellfc);
426 setup_bonded_threading(fr, &top->idef);
428 if (ir->pull && PAR(cr))
430 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
434 if (DOMAINDECOMP(cr))
436 /* Distribute the charge groups over the nodes from the master node */
437 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
438 state_global, top_global, ir,
439 state, &f, mdatoms, top, fr,
440 vsite, shellfc, constr,
441 nrnb, wcycle, FALSE);
445 update_mdatoms(mdatoms, state->lambda[efptMASS]);
447 if (opt2bSet("-cpi", nfile, fnm))
449 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
453 bStateFromCP = FALSE;
458 init_expanded_ensemble(bStateFromCP,ir,&mcrng,&state->dfhist);
465 /* Update mdebin with energy history if appending to output files */
466 if (Flags & MD_APPENDFILES)
468 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
472 /* We might have read an energy history from checkpoint,
473 * free the allocated memory and reset the counts.
475 done_energyhistory(&state_global->enerhist);
476 init_energyhistory(&state_global->enerhist);
479 /* Set the initial energy history in state by updating once */
480 update_energyhistory(&state_global->enerhist, mdebin);
483 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
485 /* Set the random state if we read a checkpoint file */
486 set_stochd_state(upd, state);
489 if (state->flags & (1<<estMC_RNG))
491 set_mc_state(mcrng, state);
494 /* Initialize constraints */
497 if (!DOMAINDECOMP(cr))
499 set_constraints(constr, top, ir, mdatoms, cr);
503 /* Check whether we have to GCT stuff */
504 bTCR = ftp2bSet(efGCT, nfile, fnm);
509 fprintf(stderr, "Will do General Coupling Theory!\n");
511 gnx = top_global->mols.nr;
513 for (i = 0; (i < gnx); i++)
521 /* We need to be sure replica exchange can only occur
522 * when the energies are current */
523 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
524 "repl_ex_nst", &repl_ex_nst);
525 /* This check needs to happen before inter-simulation
526 * signals are initialized, too */
528 if (repl_ex_nst > 0 && MASTER(cr))
530 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
531 repl_ex_nst, repl_ex_nex, repl_ex_seed);
534 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
535 * With perturbed charges with soft-core we should not change the cut-off.
537 if ((Flags & MD_TUNEPME) &&
538 EEL_PME(fr->eeltype) &&
539 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
540 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
543 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
545 if (cr->duty & DUTY_PME)
547 /* Start tuning right away, as we can't measure the load */
548 bPMETuneRunning = TRUE;
552 /* Separate PME nodes, we can measure the PP/PME load balance */
557 if (!ir->bContinuation && !bRerunMD)
559 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
561 /* Set the velocities of frozen particles to zero */
562 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
564 for (m = 0; m < DIM; m++)
566 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
576 /* Constrain the initial coordinates and velocities */
577 do_constrain_first(fplog, constr, ir, mdatoms, state, f,
578 graph, cr, nrnb, fr, top, shake_vir);
582 /* Construct the virtual sites for the initial configuration */
583 construct_vsites(fplog, vsite, state->x, nrnb, ir->delta_t, NULL,
584 top->idef.iparams, top->idef.il,
585 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
591 /* set free energy calculation frequency as the minimum
592 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
593 nstfep = ir->fepvals->nstdhdl;
596 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl,nstfep);
600 nstfep = gmx_greatest_common_divisor(repl_ex_nst,nstfep);
603 /* I'm assuming we need global communication the first time! MRS */
604 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
605 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
606 | (bVV ? CGLO_PRESSURE : 0)
607 | (bVV ? CGLO_CONSTRAINT : 0)
608 | (bRerunMD ? CGLO_RERUNMD : 0)
609 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
611 bSumEkinhOld = FALSE;
612 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
613 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
614 constr, NULL, FALSE, state->box,
615 top_global, &pcurr, top_global->natoms, &bSumEkinhOld, cglo_flags);
616 if (ir->eI == eiVVAK)
618 /* a second call to get the half step temperature initialized as well */
619 /* we do the same call as above, but turn the pressure off -- internally to
620 compute_globals, this is recognized as a velocity verlet half-step
621 kinetic energy calculation. This minimized excess variables, but
622 perhaps loses some logic?*/
624 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
625 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
626 constr, NULL, FALSE, state->box,
627 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
628 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
631 /* Calculate the initial half step temperature, and save the ekinh_old */
632 if (!(Flags & MD_STARTFROMCPT))
634 for (i = 0; (i < ir->opts.ngtc); i++)
636 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
641 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
642 and there is no previous step */
645 /* if using an iterative algorithm, we need to create a working directory for the state. */
648 bufstate = init_bufstate(state);
652 snew(xcopy, state->natoms);
653 snew(vcopy, state->natoms);
654 copy_rvecn(state->x, xcopy, 0, state->natoms);
655 copy_rvecn(state->v, vcopy, 0, state->natoms);
656 copy_mat(state->box, boxcopy);
659 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
660 temperature control */
661 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
665 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
668 "RMS relative constraint deviation after constraining: %.2e\n",
669 constr_rmsd(constr, FALSE));
671 if (EI_STATE_VELOCITY(ir->eI))
673 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
677 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
678 " input trajectory '%s'\n\n",
679 *(top_global->name), opt2fn("-rerun", nfile, fnm));
682 fprintf(stderr, "Calculated time to finish depends on nsteps from "
683 "run input file,\nwhich may not correspond to the time "
684 "needed to process input trajectory.\n\n");
690 fprintf(stderr, "starting mdrun '%s'\n",
691 *(top_global->name));
694 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
698 sprintf(tbuf, "%s", "infinite");
700 if (ir->init_step > 0)
702 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
703 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
704 gmx_step_str(ir->init_step, sbuf2),
705 ir->init_step*ir->delta_t);
709 fprintf(stderr, "%s steps, %s ps.\n",
710 gmx_step_str(ir->nsteps, sbuf), tbuf);
713 fprintf(fplog, "\n");
716 print_start(fplog, cr, runtime, "mdrun");
717 runtime_start(runtime);
718 wallcycle_start(wcycle, ewcRUN);
720 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
722 chkpt_ret = fcCheckPointParallel( cr->nodeid,
726 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
731 /***********************************************************
735 ************************************************************/
737 /* if rerunMD then read coordinates and velocities from input trajectory */
740 if (getenv("GMX_FORCE_UPDATE"))
748 bNotLastFrame = read_first_frame(oenv, &status,
749 opt2fn("-rerun", nfile, fnm),
750 &rerun_fr, TRX_NEED_X | TRX_READ_V);
751 if (rerun_fr.natoms != top_global->natoms)
754 "Number of atoms in trajectory (%d) does not match the "
755 "run input file (%d)\n",
756 rerun_fr.natoms, top_global->natoms);
758 if (ir->ePBC != epbcNONE)
762 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);
764 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
766 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
773 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
776 if (ir->ePBC != epbcNONE)
778 /* Set the shift vectors.
779 * Necessary here when have a static box different from the tpr box.
781 calc_shifts(rerun_fr.box, fr->shift_vec);
785 /* loop over MD steps or if rerunMD to end of input trajectory */
787 /* Skip the first Nose-Hoover integration when we get the state from tpx */
788 bStateFromTPX = !bStateFromCP;
789 bInitStep = bFirstStep && (bStateFromTPX || bVV);
790 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
792 bSumEkinhOld = FALSE;
795 init_global_signals(&gs, cr, ir, repl_ex_nst);
797 step = ir->init_step;
800 if (ir->nstlist == -1)
802 init_nlistheuristics(&nlh, bGStatEveryStep, step);
805 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
807 /* check how many steps are left in other sims */
808 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
812 /* and stop now if we should */
813 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
814 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
815 while (!bLastStep || (bRerunMD && bNotLastFrame))
818 wallcycle_start(wcycle, ewcSTEP);
820 GMX_MPE_LOG(ev_timestep1);
826 step = rerun_fr.step;
827 step_rel = step - ir->init_step;
840 bLastStep = (step_rel == ir->nsteps);
841 t = t0 + step*ir->delta_t;
844 if (ir->efep != efepNO || ir->bSimTemp)
846 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
847 requiring different logic. */
849 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
850 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
851 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
852 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
853 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
858 update_annealing_target_temp(&(ir->opts), t);
863 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
865 for (i = 0; i < state_global->natoms; i++)
867 copy_rvec(rerun_fr.x[i], state_global->x[i]);
871 for (i = 0; i < state_global->natoms; i++)
873 copy_rvec(rerun_fr.v[i], state_global->v[i]);
878 for (i = 0; i < state_global->natoms; i++)
880 clear_rvec(state_global->v[i]);
884 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
885 " Ekin, temperature and pressure are incorrect,\n"
886 " the virial will be incorrect when constraints are present.\n"
888 bRerunWarnNoV = FALSE;
892 copy_mat(rerun_fr.box, state_global->box);
893 copy_mat(state_global->box, state->box);
895 if (vsite && (Flags & MD_RERUN_VSITE))
897 if (DOMAINDECOMP(cr))
899 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
903 /* Following is necessary because the graph may get out of sync
904 * with the coordinates if we only have every N'th coordinate set
906 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
907 shift_self(graph, state->box, state->x);
909 construct_vsites(fplog, vsite, state->x, nrnb, ir->delta_t, state->v,
910 top->idef.iparams, top->idef.il,
911 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
914 unshift_self(graph, state->box, state->x);
919 /* Stop Center of Mass motion */
920 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
922 /* Copy back starting coordinates in case we're doing a forcefield scan */
925 for (ii = 0; (ii < state->natoms); ii++)
927 copy_rvec(xcopy[ii], state->x[ii]);
928 copy_rvec(vcopy[ii], state->v[ii]);
930 copy_mat(boxcopy, state->box);
935 /* for rerun MD always do Neighbour Searching */
936 bNS = (bFirstStep || ir->nstlist != 0);
941 /* Determine whether or not to do Neighbour Searching and LR */
942 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
944 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
945 (ir->nstlist == -1 && nlh.nabnsb > 0));
947 if (bNS && ir->nstlist == -1)
949 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
953 /* check whether we should stop because another simulation has
957 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
958 (multisim_nsteps != ir->nsteps) )
965 "Stopping simulation %d because another one has finished\n",
969 gs.sig[eglsCHKPT] = 1;
974 /* < 0 means stop at next step, > 0 means stop at next NS step */
975 if ( (gs.set[eglsSTOPCOND] < 0) ||
976 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
981 /* Determine whether or not to update the Born radii if doing GB */
982 bBornRadii = bFirstStep;
983 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
988 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
989 do_verbose = bVerbose &&
990 (step % stepout == 0 || bFirstStep || bLastStep);
992 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
1000 bMasterState = FALSE;
1001 /* Correct the new box if it is too skewed */
1002 if (DYNAMIC_BOX(*ir))
1004 if (correct_box(fplog, step, state->box, graph))
1006 bMasterState = TRUE;
1009 if (DOMAINDECOMP(cr) && bMasterState)
1011 dd_collect_state(cr->dd, state, state_global);
1015 if (DOMAINDECOMP(cr))
1017 /* Repartition the domain decomposition */
1018 wallcycle_start(wcycle, ewcDOMDEC);
1019 dd_partition_system(fplog, step, cr,
1020 bMasterState, nstglobalcomm,
1021 state_global, top_global, ir,
1022 state, &f, mdatoms, top, fr,
1023 vsite, shellfc, constr,
1025 do_verbose && !bPMETuneRunning);
1026 wallcycle_stop(wcycle, ewcDOMDEC);
1027 /* If using an iterative integrator, reallocate space to match the decomposition */
1031 if (MASTER(cr) && do_log && !bFFscan)
1033 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
1036 if (ir->efep != efepNO)
1038 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1041 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1044 /* We need the kinetic energy at minus the half step for determining
1045 * the full step kinetic energy and possibly for T-coupling.*/
1046 /* This may not be quite working correctly yet . . . . */
1047 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1048 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1049 constr, NULL, FALSE, state->box,
1050 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1051 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1053 clear_mat(force_vir);
1055 /* Ionize the atoms if necessary */
1058 ionize(fplog, oenv, mdatoms, top_global, t, ir, state->x, state->v,
1059 mdatoms->start, mdatoms->start+mdatoms->homenr, state->box, cr);
1062 /* Update force field in ffscan program */
1065 if (update_forcefield(fplog,
1067 mdatoms->nr, state->x, state->box))
1075 GMX_MPE_LOG(ev_timestep2);
1077 /* We write a checkpoint at this MD step when:
1078 * either at an NS step when we signalled through gs,
1079 * or at the last step (but not when we do not want confout),
1080 * but never at the first step or with rerun.
1082 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1083 (bLastStep && (Flags & MD_CONFOUT))) &&
1084 step > ir->init_step && !bRerunMD);
1087 gs.set[eglsCHKPT] = 0;
1090 /* Determine the energy and pressure:
1091 * at nstcalcenergy steps and at energy output steps (set below).
1093 if (EI_VV(ir->eI) && (!bInitStep))
1095 /* for vv, the first half of the integration actually corresponds
1096 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1097 but the virial needs to be calculated on both the current step and the 'next' step. Future
1098 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1100 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1101 bCalcVir = bCalcEner ||
1102 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1106 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1107 bCalcVir = bCalcEner ||
1108 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1111 /* Do we need global communication ? */
1112 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1113 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1114 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1116 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1118 if (do_ene || do_log)
1125 /* these CGLO_ options remain the same throughout the iteration */
1126 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1127 (bGStat ? CGLO_GSTAT : 0)
1130 force_flags = (GMX_FORCE_STATECHANGED |
1131 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1132 GMX_FORCE_ALLFORCES |
1134 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1135 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1136 (bDoFEP ? GMX_FORCE_DHDL : 0)
1141 if (do_per_step(step, ir->nstcalclr))
1143 force_flags |= GMX_FORCE_DO_LR;
1149 /* Now is the time to relax the shells */
1150 count = relax_shell_flexcon(fplog, cr, bVerbose, bFFscan ? step+1 : step,
1151 ir, bNS, force_flags,
1152 bStopCM, top, top_global,
1154 state, f, force_vir, mdatoms,
1155 nrnb, wcycle, graph, groups,
1156 shellfc, fr, bBornRadii, t, mu_tot,
1157 state->natoms, &bConverged, vsite,
1168 /* The coordinates (x) are shifted (to get whole molecules)
1170 * This is parallellized as well, and does communication too.
1171 * Check comments in sim_util.c
1173 do_force(fplog, cr, ir, step, nrnb, wcycle, top, top_global, groups,
1174 state->box, state->x, &state->hist,
1175 f, force_vir, mdatoms, enerd, fcd,
1176 state->lambda, graph,
1177 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1178 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1181 GMX_BARRIER(cr->mpi_comm_mygroup);
1185 mu_aver = calc_mu_aver(cr, state->x, mdatoms->chargeA,
1186 mu_tot, &top_global->mols, mdatoms, gnx, grpindex);
1189 if (bTCR && bFirstStep)
1191 tcr = init_coupling(fplog, nfile, fnm, cr, fr, mdatoms, &(top->idef));
1192 fprintf(fplog, "Done init_coupling\n");
1196 if (bVV && !bStartingFromCpt && !bRerunMD)
1197 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1199 if (ir->eI == eiVV && bInitStep)
1201 /* if using velocity verlet with full time step Ekin,
1202 * take the first half step only to compute the
1203 * virial for the first step. From there,
1204 * revert back to the initial coordinates
1205 * so that the input is actually the initial step.
1207 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1211 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1212 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1215 /* If we are using twin-range interactions where the long-range component
1216 * is only evaluated every nstcalclr>1 steps, we should do a special update
1217 * step to combine the long-range forces on these steps.
1218 * For nstcalclr=1 this is not done, since the forces would have been added
1219 * directly to the short-range forces already.
1221 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1223 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1224 f, bUpdateDoLR, fr->f_twin, fcd,
1225 ekind, M, wcycle, upd, bInitStep, etrtVELOCITY1,
1226 cr, nrnb, constr, &top->idef);
1228 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1230 gmx_iterate_init(&iterate, TRUE);
1232 /* for iterations, we save these vectors, as we will be self-consistently iterating
1235 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1237 /* save the state */
1238 if (iterate.bIterationActive)
1240 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1243 bFirstIterate = TRUE;
1244 while (bFirstIterate || iterate.bIterationActive)
1246 if (iterate.bIterationActive)
1248 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1249 if (bFirstIterate && bTrotter)
1251 /* The first time through, we need a decent first estimate
1252 of veta(t+dt) to compute the constraints. Do
1253 this by computing the box volume part of the
1254 trotter integration at this time. Nothing else
1255 should be changed by this routine here. If
1256 !(first time), we start with the previous value
1259 veta_save = state->veta;
1260 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1261 vetanew = state->veta;
1262 state->veta = veta_save;
1267 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1269 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1270 state, fr->bMolPBC, graph, f,
1271 &top->idef, shake_vir, NULL,
1272 cr, nrnb, wcycle, upd, constr,
1273 bInitStep, TRUE, bCalcVir, vetanew);
1275 if (!bOK && !bFFscan)
1277 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1283 /* Need to unshift here if a do_force has been
1284 called in the previous step */
1285 unshift_self(graph, state->box, state->x);
1288 /* if VV, compute the pressure and constraints */
1289 /* For VV2, we strictly only need this if using pressure
1290 * control, but we really would like to have accurate pressures
1292 * Think about ways around this in the future?
1293 * For now, keep this choice in comments.
1295 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1296 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1298 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1299 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1301 bSumEkinhOld = TRUE;
1303 /* for vv, the first half of the integration actually corresponds to the previous step.
1304 So we need information from the last step in the first half of the integration */
1305 if (bGStat || do_per_step(step-1, nstglobalcomm))
1307 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1308 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1309 constr, NULL, FALSE, state->box,
1310 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1313 | (bTemp ? CGLO_TEMPERATURE : 0)
1314 | (bPres ? CGLO_PRESSURE : 0)
1315 | (bPres ? CGLO_CONSTRAINT : 0)
1316 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1317 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1320 /* explanation of above:
1321 a) We compute Ekin at the full time step
1322 if 1) we are using the AveVel Ekin, and it's not the
1323 initial step, or 2) if we are using AveEkin, but need the full
1324 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1325 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1326 EkinAveVel because it's needed for the pressure */
1328 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1333 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1334 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1341 /* We need the kinetic energy at minus the half step for determining
1342 * the full step kinetic energy and possibly for T-coupling.*/
1343 /* This may not be quite working correctly yet . . . . */
1344 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1345 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1346 constr, NULL, FALSE, state->box,
1347 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1348 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1353 if (iterate.bIterationActive &&
1354 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1355 state->veta, &vetanew))
1359 bFirstIterate = FALSE;
1362 if (bTrotter && !bInitStep)
1364 copy_mat(shake_vir, state->svir_prev);
1365 copy_mat(force_vir, state->fvir_prev);
1366 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1368 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1369 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE, FALSE);
1370 enerd->term[F_EKIN] = trace(ekind->ekin);
1373 /* if it's the initial step, we performed this first step just to get the constraint virial */
1374 if (bInitStep && ir->eI == eiVV)
1376 copy_rvecn(cbuf, state->v, 0, state->natoms);
1379 GMX_MPE_LOG(ev_timestep1);
1382 /* MRS -- now done iterating -- compute the conserved quantity */
1385 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1388 last_ekin = enerd->term[F_EKIN];
1390 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1392 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1394 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1397 sum_dhdl(enerd, state->lambda, ir->fepvals);
1401 /* ######## END FIRST UPDATE STEP ############## */
1402 /* ######## If doing VV, we now have v(dt) ###### */
1405 /* perform extended ensemble sampling in lambda - we don't
1406 actually move to the new state before outputting
1407 statistics, but if performing simulated tempering, we
1408 do update the velocities and the tau_t. */
1410 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, mcrng, state->v, mdatoms);
1411 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1412 copy_df_history(&state_global->dfhist,&state->dfhist);
1414 /* ################## START TRAJECTORY OUTPUT ################# */
1416 /* Now we have the energies and forces corresponding to the
1417 * coordinates at time t. We must output all of this before
1419 * for RerunMD t is read from input trajectory
1421 GMX_MPE_LOG(ev_output_start);
1424 if (do_per_step(step, ir->nstxout))
1426 mdof_flags |= MDOF_X;
1428 if (do_per_step(step, ir->nstvout))
1430 mdof_flags |= MDOF_V;
1432 if (do_per_step(step, ir->nstfout))
1434 mdof_flags |= MDOF_F;
1436 if (do_per_step(step, ir->nstxtcout))
1438 mdof_flags |= MDOF_XTC;
1442 mdof_flags |= MDOF_CPT;
1446 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1449 /* Enforce writing positions and velocities at end of run */
1450 mdof_flags |= (MDOF_X | MDOF_V);
1456 fcReportProgress( ir->nsteps, step );
1459 #if defined(__native_client__)
1460 fcCheckin(MASTER(cr));
1463 /* sync bCPT and fc record-keeping */
1464 if (bCPT && MASTER(cr))
1466 fcRequestCheckPoint();
1470 if (mdof_flags != 0)
1472 wallcycle_start(wcycle, ewcTRAJ);
1475 if (state->flags & (1<<estLD_RNG))
1477 get_stochd_state(upd, state);
1479 if (state->flags & (1<<estMC_RNG))
1481 get_mc_state(mcrng, state);
1487 state_global->ekinstate.bUpToDate = FALSE;
1491 update_ekinstate(&state_global->ekinstate, ekind);
1492 state_global->ekinstate.bUpToDate = TRUE;
1494 update_energyhistory(&state_global->enerhist, mdebin);
1497 write_traj(fplog, cr, outf, mdof_flags, top_global,
1498 step, t, state, state_global, f, f_global, &n_xtc, &x_xtc);
1505 if (bLastStep && step_rel == ir->nsteps &&
1506 (Flags & MD_CONFOUT) && MASTER(cr) &&
1507 !bRerunMD && !bFFscan)
1509 /* x and v have been collected in write_traj,
1510 * because a checkpoint file will always be written
1513 fprintf(stderr, "\nWriting final coordinates.\n");
1516 /* Make molecules whole only for confout writing */
1517 do_pbc_mtop(fplog, ir->ePBC, state->box, top_global, state_global->x);
1519 write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm),
1520 *top_global->name, top_global,
1521 state_global->x, state_global->v,
1522 ir->ePBC, state->box);
1525 wallcycle_stop(wcycle, ewcTRAJ);
1527 GMX_MPE_LOG(ev_output_finish);
1529 /* kludge -- virial is lost with restart for NPT control. Must restart */
1530 if (bStartingFromCpt && bVV)
1532 copy_mat(state->svir_prev, shake_vir);
1533 copy_mat(state->fvir_prev, force_vir);
1535 /* ################## END TRAJECTORY OUTPUT ################ */
1537 /* Determine the wallclock run time up till now */
1538 run_time = gmx_gettime() - (double)runtime->real;
1540 /* Check whether everything is still allright */
1541 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1542 #ifdef GMX_THREAD_MPI
1547 /* this is just make gs.sig compatible with the hack
1548 of sending signals around by MPI_Reduce with together with
1550 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1552 gs.sig[eglsSTOPCOND] = 1;
1554 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1556 gs.sig[eglsSTOPCOND] = -1;
1558 /* < 0 means stop at next step, > 0 means stop at next NS step */
1562 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1563 gmx_get_signal_name(),
1564 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1568 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1569 gmx_get_signal_name(),
1570 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1572 handled_stop_condition = (int)gmx_get_stop_condition();
1574 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1575 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1576 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1578 /* Signal to terminate the run */
1579 gs.sig[eglsSTOPCOND] = 1;
1582 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1584 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1587 if (bResetCountersHalfMaxH && MASTER(cr) &&
1588 run_time > max_hours*60.0*60.0*0.495)
1590 gs.sig[eglsRESETCOUNTERS] = 1;
1593 if (ir->nstlist == -1 && !bRerunMD)
1595 /* When bGStatEveryStep=FALSE, global_stat is only called
1596 * when we check the atom displacements, not at NS steps.
1597 * This means that also the bonded interaction count check is not
1598 * performed immediately after NS. Therefore a few MD steps could
1599 * be performed with missing interactions.
1600 * But wrong energies are never written to file,
1601 * since energies are only written after global_stat
1604 if (step >= nlh.step_nscheck)
1606 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1607 nlh.scale_tot, state->x);
1611 /* This is not necessarily true,
1612 * but step_nscheck is determined quite conservatively.
1618 /* In parallel we only have to check for checkpointing in steps
1619 * where we do global communication,
1620 * otherwise the other nodes don't know.
1622 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1625 run_time >= nchkpt*cpt_period*60.0)) &&
1626 gs.set[eglsCHKPT] == 0)
1628 gs.sig[eglsCHKPT] = 1;
1631 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1636 update_tcouple(fplog, step, ir, state, ekind, wcycle, upd, &MassQ, mdatoms);
1638 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1640 gmx_bool bIfRandomize;
1641 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr, DOMAINDECOMP(cr));
1642 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1643 if (constr && bIfRandomize)
1645 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1646 state, fr->bMolPBC, graph, f,
1647 &top->idef, tmp_vir, NULL,
1648 cr, nrnb, wcycle, upd, constr,
1649 bInitStep, TRUE, bCalcVir, vetanew);
1654 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1656 gmx_iterate_init(&iterate, TRUE);
1657 /* for iterations, we save these vectors, as we will be redoing the calculations */
1658 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1661 bFirstIterate = TRUE;
1662 while (bFirstIterate || iterate.bIterationActive)
1664 /* We now restore these vectors to redo the calculation with improved extended variables */
1665 if (iterate.bIterationActive)
1667 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1670 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1671 so scroll down for that logic */
1673 /* ######### START SECOND UPDATE STEP ################# */
1674 GMX_MPE_LOG(ev_update_start);
1675 /* Box is changed in update() when we do pressure coupling,
1676 * but we should still use the old box for energy corrections and when
1677 * writing it to the energy file, so it matches the trajectory files for
1678 * the same timestep above. Make a copy in a separate array.
1680 copy_mat(state->box, lastbox);
1685 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1687 wallcycle_start(wcycle, ewcUPDATE);
1688 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1691 if (iterate.bIterationActive)
1699 /* we use a new value of scalevir to converge the iterations faster */
1700 scalevir = tracevir/trace(shake_vir);
1702 msmul(shake_vir, scalevir, shake_vir);
1703 m_add(force_vir, shake_vir, total_vir);
1704 clear_mat(shake_vir);
1706 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1707 /* We can only do Berendsen coupling after we have summed
1708 * the kinetic energy or virial. Since the happens
1709 * in global_state after update, we should only do it at
1710 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1715 update_tcouple(fplog, step, ir, state, ekind, wcycle, upd, &MassQ, mdatoms);
1716 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, wcycle,
1722 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1724 /* velocity half-step update */
1725 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1726 bUpdateDoLR, fr->f_twin, fcd,
1727 ekind, M, wcycle, upd, FALSE, etrtVELOCITY2,
1728 cr, nrnb, constr, &top->idef);
1731 /* Above, initialize just copies ekinh into ekin,
1732 * it doesn't copy position (for VV),
1733 * and entire integrator for MD.
1736 if (ir->eI == eiVVAK)
1738 copy_rvecn(state->x, cbuf, 0, state->natoms);
1740 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1742 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1743 bUpdateDoLR, fr->f_twin, fcd,
1744 ekind, M, wcycle, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1745 wallcycle_stop(wcycle, ewcUPDATE);
1747 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1748 fr->bMolPBC, graph, f,
1749 &top->idef, shake_vir, force_vir,
1750 cr, nrnb, wcycle, upd, constr,
1751 bInitStep, FALSE, bCalcVir, state->veta);
1753 if (ir->eI == eiVVAK)
1755 /* erase F_EKIN and F_TEMP here? */
1756 /* just compute the kinetic energy at the half step to perform a trotter step */
1757 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1758 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1759 constr, NULL, FALSE, lastbox,
1760 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1761 cglo_flags | CGLO_TEMPERATURE
1763 wallcycle_start(wcycle, ewcUPDATE);
1764 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1765 /* now we know the scaling, we can compute the positions again again */
1766 copy_rvecn(cbuf, state->x, 0, state->natoms);
1768 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1770 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1771 bUpdateDoLR, fr->f_twin, fcd,
1772 ekind, M, wcycle, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1773 wallcycle_stop(wcycle, ewcUPDATE);
1775 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1776 /* are the small terms in the shake_vir here due
1777 * to numerical errors, or are they important
1778 * physically? I'm thinking they are just errors, but not completely sure.
1779 * For now, will call without actually constraining, constr=NULL*/
1780 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1781 state, fr->bMolPBC, graph, f,
1782 &top->idef, tmp_vir, force_vir,
1783 cr, nrnb, wcycle, upd, NULL,
1784 bInitStep, FALSE, bCalcVir,
1787 if (!bOK && !bFFscan)
1789 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1792 if (fr->bSepDVDL && fplog && do_log)
1794 fprintf(fplog, sepdvdlformat, "Constraint dV/dl", 0.0, dvdl_constr);
1798 /* this factor or 2 correction is necessary
1799 because half of the constraint force is removed
1800 in the vv step, so we have to double it. See
1801 the Redmine issue #1255. It is not yet clear
1802 if the factor of 2 is exact, or just a very
1803 good approximation, and this will be
1804 investigated. The next step is to see if this
1805 can be done adding a dhdl contribution from the
1806 rattle step, but this is somewhat more
1807 complicated with the current code. Will be
1808 investigated, hopefully for 4.6.3. However,
1809 this current solution is much better than
1810 having it completely wrong.
1812 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1816 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1821 /* Need to unshift here */
1822 unshift_self(graph, state->box, state->x);
1825 GMX_BARRIER(cr->mpi_comm_mygroup);
1826 GMX_MPE_LOG(ev_update_finish);
1830 wallcycle_start(wcycle, ewcVSITECONSTR);
1833 shift_self(graph, state->box, state->x);
1835 construct_vsites(fplog, vsite, state->x, nrnb, ir->delta_t, state->v,
1836 top->idef.iparams, top->idef.il,
1837 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1841 unshift_self(graph, state->box, state->x);
1843 wallcycle_stop(wcycle, ewcVSITECONSTR);
1846 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1847 /* With Leap-Frog we can skip compute_globals at
1848 * non-communication steps, but we need to calculate
1849 * the kinetic energy one step before communication.
1851 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1853 if (ir->nstlist == -1 && bFirstIterate)
1855 gs.sig[eglsNABNSB] = nlh.nabnsb;
1857 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1858 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1860 bFirstIterate ? &gs : NULL,
1861 (step_rel % gs.nstms == 0) &&
1862 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1864 top_global, &pcurr, top_global->natoms, &bSumEkinhOld,
1866 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1867 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1868 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1869 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1870 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1871 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1874 if (ir->nstlist == -1 && bFirstIterate)
1876 nlh.nabnsb = gs.set[eglsNABNSB];
1877 gs.set[eglsNABNSB] = 0;
1880 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1881 /* ############# END CALC EKIN AND PRESSURE ################# */
1883 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1884 the virial that should probably be addressed eventually. state->veta has better properies,
1885 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1886 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1888 if (iterate.bIterationActive &&
1889 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1890 trace(shake_vir), &tracevir))
1894 bFirstIterate = FALSE;
1897 if (!bVV || bRerunMD)
1899 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1900 sum_dhdl(enerd, state->lambda, ir->fepvals);
1902 update_box(fplog, step, ir, mdatoms, state, graph, f,
1903 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, wcycle, upd, bInitStep, FALSE);
1905 /* ################# END UPDATE STEP 2 ################# */
1906 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1908 /* The coordinates (x) were unshifted in update */
1909 if (bFFscan && (shellfc == NULL || bConverged))
1911 if (print_forcefield(fplog, enerd->term, mdatoms->homenr,
1913 &(top_global->mols), mdatoms->massT, pres))
1917 fprintf(stderr, "\n");
1923 /* We will not sum ekinh_old,
1924 * so signal that we still have to do it.
1926 bSumEkinhOld = TRUE;
1931 /* Only do GCT when the relaxation of shells (minimization) has converged,
1932 * otherwise we might be coupling to bogus energies.
1933 * In parallel we must always do this, because the other sims might
1937 /* Since this is called with the new coordinates state->x, I assume
1938 * we want the new box state->box too. / EL 20040121
1940 do_coupling(fplog, oenv, nfile, fnm, tcr, t, step, enerd->term, fr,
1942 mdatoms, &(top->idef), mu_aver,
1943 top_global->mols.nr, cr,
1944 state->box, total_vir, pres,
1945 mu_tot, state->x, f, bConverged);
1949 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1951 /* use the directly determined last velocity, not actually the averaged half steps */
1952 if (bTrotter && ir->eI == eiVV)
1954 enerd->term[F_EKIN] = last_ekin;
1956 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1960 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1964 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1966 /* Check for excessively large energies */
1970 real etot_max = 1e200;
1972 real etot_max = 1e30;
1974 if (fabs(enerd->term[F_ETOT]) > etot_max)
1976 fprintf(stderr, "Energy too large (%g), giving up\n",
1977 enerd->term[F_ETOT]);
1980 /* ######### END PREPARING EDR OUTPUT ########### */
1982 /* Time for performance */
1983 if (((step % stepout) == 0) || bLastStep)
1985 runtime_upd_proc(runtime);
1991 gmx_bool do_dr, do_or;
1993 if (fplog && do_log && bDoExpanded)
1995 /* only needed if doing expanded ensemble */
1996 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1997 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1999 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
2003 upd_mdebin(mdebin, bDoDHDL, TRUE,
2004 t, mdatoms->tmass, enerd, state,
2005 ir->fepvals, ir->expandedvals, lastbox,
2006 shake_vir, force_vir, total_vir, pres,
2007 ekind, mu_tot, constr);
2011 upd_mdebin_step(mdebin);
2014 do_dr = do_per_step(step, ir->nstdisreout);
2015 do_or = do_per_step(step, ir->nstorireout);
2017 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
2019 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
2021 if (ir->ePull != epullNO)
2023 pull_print_output(ir->pull, step, t);
2026 if (do_per_step(step, ir->nstlog))
2028 if (fflush(fplog) != 0)
2030 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
2036 /* Have to do this part _after_ outputting the logfile and the edr file */
2037 /* Gets written into the state at the beginning of next loop*/
2038 state->fep_state = lamnew;
2041 /* Remaining runtime */
2042 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
2046 fprintf(stderr, "\n");
2048 print_time(stderr, runtime, step, ir, cr);
2051 /* Replica exchange */
2053 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
2054 do_per_step(step, repl_ex_nst))
2056 bExchanged = replica_exchange(fplog, cr, repl_ex,
2057 state_global, enerd,
2060 if (bExchanged && DOMAINDECOMP(cr))
2062 dd_partition_system(fplog, step, cr, TRUE, 1,
2063 state_global, top_global, ir,
2064 state, &f, mdatoms, top, fr,
2065 vsite, shellfc, constr,
2066 nrnb, wcycle, FALSE);
2072 bStartingFromCpt = FALSE;
2074 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
2075 /* With all integrators, except VV, we need to retain the pressure
2076 * at the current step for coupling at the next step.
2078 if ((state->flags & (1<<estPRES_PREV)) &&
2080 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
2082 /* Store the pressure in t_state for pressure coupling
2083 * at the next MD step.
2085 copy_mat(pres, state->pres_prev);
2088 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2090 if ( (membed != NULL) && (!bLastStep) )
2092 rescale_membed(step_rel, membed, state_global->x);
2099 /* read next frame from input trajectory */
2100 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
2105 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
2109 if (!bRerunMD || !rerun_fr.bStep)
2111 /* increase the MD step number */
2116 cycles = wallcycle_stop(wcycle, ewcSTEP);
2117 if (DOMAINDECOMP(cr) && wcycle)
2119 dd_cycles_add(cr->dd, cycles, ddCyclStep);
2122 if (bPMETuneRunning || bPMETuneTry)
2124 /* PME grid + cut-off optimization with GPUs or PME nodes */
2126 /* Count the total cycles over the last steps */
2127 cycles_pmes += cycles;
2129 /* We can only switch cut-off at NS steps */
2130 if (step % ir->nstlist == 0)
2132 /* PME grid + cut-off optimization with GPUs or PME nodes */
2135 if (DDMASTER(cr->dd))
2137 /* PME node load is too high, start tuning */
2138 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2140 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
2142 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2144 bPMETuneTry = FALSE;
2147 if (bPMETuneRunning)
2149 /* init_step might not be a multiple of nstlist,
2150 * but the first cycle is always skipped anyhow.
2153 pme_load_balance(pme_loadbal, cr,
2154 (bVerbose && MASTER(cr)) ? stderr : NULL,
2156 ir, state, cycles_pmes,
2157 fr->ic, fr->nbv, &fr->pmedata,
2160 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2161 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2162 fr->rlist = fr->ic->rlist;
2163 fr->rlistlong = fr->ic->rlistlong;
2164 fr->rcoulomb = fr->ic->rcoulomb;
2165 fr->rvdw = fr->ic->rvdw;
2171 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2172 gs.set[eglsRESETCOUNTERS] != 0)
2174 /* Reset all the counters related to performance over the run */
2175 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, runtime,
2176 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2177 wcycle_set_reset_counters(wcycle, -1);
2178 if (!(cr->duty & DUTY_PME))
2180 /* Tell our PME node to reset its counters */
2181 gmx_pme_send_resetcounters(cr, step);
2183 /* Correct max_hours for the elapsed time */
2184 max_hours -= run_time/(60.0*60.0);
2185 bResetCountersHalfMaxH = FALSE;
2186 gs.set[eglsRESETCOUNTERS] = 0;
2190 /* End of main MD loop */
2194 runtime_end(runtime);
2196 if (bRerunMD && MASTER(cr))
2201 if (!(cr->duty & DUTY_PME))
2203 /* Tell the PME only node to finish */
2204 gmx_pme_send_finish(cr);
2209 if (ir->nstcalcenergy > 0 && !bRerunMD)
2211 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
2212 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
2220 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2222 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)));
2223 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
2226 if (pme_loadbal != NULL)
2228 pme_loadbal_done(pme_loadbal, cr, fplog,
2229 fr->nbv != NULL && fr->nbv->bUseGPU);
2232 if (shellfc && fplog)
2234 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
2235 (nconverged*100.0)/step_rel);
2236 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
2240 if (repl_ex_nst > 0 && MASTER(cr))
2242 print_replica_exchange_statistics(fplog, repl_ex);
2245 runtime->nsteps_done = step_rel;