3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
52 #include "gromacs/fileio/trnio.h"
53 #include "gromacs/fileio/xtcio.h"
55 #include "md_support.h"
56 #include "md_logging.h"
57 #include "gromacs/fileio/confio.h"
71 #include "domdec_network.h"
77 #include "compute_io.h"
79 #include "checkpoint.h"
80 #include "mtop_util.h"
81 #include "sighandler.h"
84 #include "pme_loadbal.h"
87 #include "types/nlistheuristics.h"
88 #include "types/iteratedconstraints.h"
89 #include "nbnxn_cuda_data_mgmt.h"
91 #include "gromacs/utility/gmxmpi.h"
92 #include "gromacs/fileio/trajectory_writing.h"
93 #include "gromacs/fileio/trxio.h"
94 #include "gromacs/timing/walltime_accounting.h"
100 static void reset_all_counters(FILE *fplog, t_commrec *cr,
101 gmx_large_int_t step,
102 gmx_large_int_t *step_rel, t_inputrec *ir,
103 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
104 gmx_walltime_accounting_t walltime_accounting,
105 nbnxn_cuda_ptr_t cu_nbv)
107 char sbuf[STEPSTRSIZE];
109 /* Reset all the counters related to performance over the run */
110 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
111 gmx_step_str(step, sbuf));
115 nbnxn_cuda_reset_timings(cu_nbv);
118 wallcycle_stop(wcycle, ewcRUN);
119 wallcycle_reset_all(wcycle);
120 if (DOMAINDECOMP(cr))
122 reset_dd_statistics_counters(cr->dd);
125 ir->init_step += *step_rel;
126 ir->nsteps -= *step_rel;
128 wallcycle_start(wcycle, ewcRUN);
129 walltime_accounting_start(walltime_accounting);
130 print_date_and_time(fplog, cr->nodeid, "Restarted time", walltime_accounting);
133 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
134 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
136 gmx_vsite_t *vsite, gmx_constr_t constr,
137 int stepout, t_inputrec *ir,
138 gmx_mtop_t *top_global,
140 t_state *state_global,
142 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
143 gmx_edsam_t ed, t_forcerec *fr,
144 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
145 real cpt_period, real max_hours,
146 const char gmx_unused *deviceOptions,
148 gmx_walltime_accounting_t walltime_accounting)
151 gmx_large_int_t step, step_rel;
153 double t, t0, lam0[efptNR];
154 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
155 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
156 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
157 bBornRadii, bStartingFromCpt;
158 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
159 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
160 bForceUpdate = FALSE, bCPT;
161 gmx_bool bMasterState;
162 int force_flags, cglo_flags;
163 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
168 t_state *bufstate = NULL;
169 matrix *scale_tot, pcoupl_mu, M, ebox;
172 gmx_repl_ex_t repl_ex = NULL;
175 t_mdebin *mdebin = NULL;
176 t_state *state = NULL;
177 rvec *f_global = NULL;
178 gmx_enerdata_t *enerd;
180 gmx_global_stat_t gstat;
181 gmx_update_t upd = NULL;
182 t_graph *graph = NULL;
184 gmx_rng_t mcrng = NULL;
185 gmx_groups_t *groups;
186 gmx_ekindata_t *ekind, *ekind_save;
187 gmx_shellfc_t shellfc;
188 int count, nconverged = 0;
191 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
193 gmx_bool bResetCountersHalfMaxH = FALSE;
194 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
195 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_large_int_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");
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,
327 top_global, n_flexible_constraints(constr),
328 (ir->bContinuation ||
329 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
330 NULL : state_global->x);
334 #ifdef GMX_THREAD_MPI
335 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
337 set_deform_reference_box(upd,
338 deform_init_init_step_tpx,
339 deform_init_box_tpx);
340 #ifdef GMX_THREAD_MPI
341 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
346 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
347 if ((io > 2000) && MASTER(cr))
350 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
355 if (DOMAINDECOMP(cr))
357 top = dd_init_local_top(top_global);
360 dd_init_local_state(cr->dd, state_global, state);
362 if (DDMASTER(cr->dd) && ir->nstfout)
364 snew(f_global, state_global->natoms);
371 /* Initialize the particle decomposition and split the topology */
372 top = split_system(fplog, top_global, ir, cr);
374 pd_cg_range(cr, &fr->cg0, &fr->hcg);
375 pd_at_range(cr, &a0, &a1);
379 top = gmx_mtop_generate_local_top(top_global, ir);
382 a1 = top_global->natoms;
385 forcerec_set_excl_load(fr, top, cr);
387 state = partdec_init_local_state(cr, state_global);
390 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
394 set_vsite_top(vsite, top, mdatoms, cr);
397 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
399 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
404 make_local_shells(cr, mdatoms, shellfc);
407 setup_bonded_threading(fr, &top->idef);
409 if (ir->pull && PAR(cr))
411 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
415 if (DOMAINDECOMP(cr))
417 /* Distribute the charge groups over the nodes from the master node */
418 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
419 state_global, top_global, ir,
420 state, &f, mdatoms, top, fr,
421 vsite, shellfc, constr,
422 nrnb, wcycle, FALSE);
426 update_mdatoms(mdatoms, state->lambda[efptMASS]);
428 if (opt2bSet("-cpi", nfile, fnm))
430 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
434 bStateFromCP = FALSE;
439 init_expanded_ensemble(bStateFromCP,ir,&mcrng,&state->dfhist);
446 /* Update mdebin with energy history if appending to output files */
447 if (Flags & MD_APPENDFILES)
449 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
453 /* We might have read an energy history from checkpoint,
454 * free the allocated memory and reset the counts.
456 done_energyhistory(&state_global->enerhist);
457 init_energyhistory(&state_global->enerhist);
460 /* Set the initial energy history in state by updating once */
461 update_energyhistory(&state_global->enerhist, mdebin);
464 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
466 /* Set the random state if we read a checkpoint file */
467 set_stochd_state(upd, state);
470 if (state->flags & (1<<estMC_RNG))
472 set_mc_state(mcrng, state);
475 /* Initialize constraints */
478 if (!DOMAINDECOMP(cr))
480 set_constraints(constr, top, ir, mdatoms, cr);
486 /* We need to be sure replica exchange can only occur
487 * when the energies are current */
488 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
489 "repl_ex_nst", &repl_ex_nst);
490 /* This check needs to happen before inter-simulation
491 * signals are initialized, too */
493 if (repl_ex_nst > 0 && MASTER(cr))
495 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
496 repl_ex_nst, repl_ex_nex, repl_ex_seed);
499 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
500 * With perturbed charges with soft-core we should not change the cut-off.
502 if ((Flags & MD_TUNEPME) &&
503 EEL_PME(fr->eeltype) &&
504 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
505 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
508 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
510 if (cr->duty & DUTY_PME)
512 /* Start tuning right away, as we can't measure the load */
513 bPMETuneRunning = TRUE;
517 /* Separate PME nodes, we can measure the PP/PME load balance */
522 if (!ir->bContinuation && !bRerunMD)
524 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
526 /* Set the velocities of frozen particles to zero */
527 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
529 for (m = 0; m < DIM; m++)
531 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
541 /* Constrain the initial coordinates and velocities */
542 do_constrain_first(fplog, constr, ir, mdatoms, state,
547 /* Construct the virtual sites for the initial configuration */
548 construct_vsites(vsite, state->x, ir->delta_t, NULL,
549 top->idef.iparams, top->idef.il,
550 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
556 /* set free energy calculation frequency as the minimum
557 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
558 nstfep = ir->fepvals->nstdhdl;
561 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
565 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
568 /* I'm assuming we need global communication the first time! MRS */
569 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
570 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
571 | (bVV ? CGLO_PRESSURE : 0)
572 | (bVV ? CGLO_CONSTRAINT : 0)
573 | (bRerunMD ? CGLO_RERUNMD : 0)
574 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
576 bSumEkinhOld = FALSE;
577 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
578 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
579 constr, NULL, FALSE, state->box,
580 top_global, &bSumEkinhOld, cglo_flags);
581 if (ir->eI == eiVVAK)
583 /* a second call to get the half step temperature initialized as well */
584 /* we do the same call as above, but turn the pressure off -- internally to
585 compute_globals, this is recognized as a velocity verlet half-step
586 kinetic energy calculation. This minimized excess variables, but
587 perhaps loses some logic?*/
589 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
590 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
591 constr, NULL, FALSE, state->box,
592 top_global, &bSumEkinhOld,
593 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
596 /* Calculate the initial half step temperature, and save the ekinh_old */
597 if (!(Flags & MD_STARTFROMCPT))
599 for (i = 0; (i < ir->opts.ngtc); i++)
601 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
606 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
607 and there is no previous step */
610 /* if using an iterative algorithm, we need to create a working directory for the state. */
613 bufstate = init_bufstate(state);
616 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
617 temperature control */
618 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
622 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
625 "RMS relative constraint deviation after constraining: %.2e\n",
626 constr_rmsd(constr, FALSE));
628 if (EI_STATE_VELOCITY(ir->eI))
630 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
634 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
635 " input trajectory '%s'\n\n",
636 *(top_global->name), opt2fn("-rerun", nfile, fnm));
639 fprintf(stderr, "Calculated time to finish depends on nsteps from "
640 "run input file,\nwhich may not correspond to the time "
641 "needed to process input trajectory.\n\n");
647 fprintf(stderr, "starting mdrun '%s'\n",
648 *(top_global->name));
651 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
655 sprintf(tbuf, "%s", "infinite");
657 if (ir->init_step > 0)
659 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
660 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
661 gmx_step_str(ir->init_step, sbuf2),
662 ir->init_step*ir->delta_t);
666 fprintf(stderr, "%s steps, %s ps.\n",
667 gmx_step_str(ir->nsteps, sbuf), tbuf);
670 fprintf(fplog, "\n");
673 /* Set and write start time */
674 walltime_accounting_start(walltime_accounting);
675 print_date_and_time(fplog, cr->nodeid, "Started mdrun", walltime_accounting);
676 wallcycle_start(wcycle, ewcRUN);
679 fprintf(fplog, "\n");
682 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
684 chkpt_ret = fcCheckPointParallel( cr->nodeid,
688 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
693 /***********************************************************
697 ************************************************************/
699 /* if rerunMD then read coordinates and velocities from input trajectory */
702 if (getenv("GMX_FORCE_UPDATE"))
710 bNotLastFrame = read_first_frame(oenv, &status,
711 opt2fn("-rerun", nfile, fnm),
712 &rerun_fr, TRX_NEED_X | TRX_READ_V);
713 if (rerun_fr.natoms != top_global->natoms)
716 "Number of atoms in trajectory (%d) does not match the "
717 "run input file (%d)\n",
718 rerun_fr.natoms, top_global->natoms);
720 if (ir->ePBC != epbcNONE)
724 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);
726 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
728 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
735 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
738 if (ir->ePBC != epbcNONE)
740 /* Set the shift vectors.
741 * Necessary here when have a static box different from the tpr box.
743 calc_shifts(rerun_fr.box, fr->shift_vec);
747 /* loop over MD steps or if rerunMD to end of input trajectory */
749 /* Skip the first Nose-Hoover integration when we get the state from tpx */
750 bStateFromTPX = !bStateFromCP;
751 bInitStep = bFirstStep && (bStateFromTPX || bVV);
752 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
754 bSumEkinhOld = FALSE;
757 init_global_signals(&gs, cr, ir, repl_ex_nst);
759 step = ir->init_step;
762 if (ir->nstlist == -1)
764 init_nlistheuristics(&nlh, bGStatEveryStep, step);
767 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
769 /* check how many steps are left in other sims */
770 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
774 /* and stop now if we should */
775 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
776 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
777 while (!bLastStep || (bRerunMD && bNotLastFrame))
780 wallcycle_start(wcycle, ewcSTEP);
786 step = rerun_fr.step;
787 step_rel = step - ir->init_step;
800 bLastStep = (step_rel == ir->nsteps);
801 t = t0 + step*ir->delta_t;
804 if (ir->efep != efepNO || ir->bSimTemp)
806 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
807 requiring different logic. */
809 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
810 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
811 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
812 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
813 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
818 update_annealing_target_temp(&(ir->opts), t);
823 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
825 for (i = 0; i < state_global->natoms; i++)
827 copy_rvec(rerun_fr.x[i], state_global->x[i]);
831 for (i = 0; i < state_global->natoms; i++)
833 copy_rvec(rerun_fr.v[i], state_global->v[i]);
838 for (i = 0; i < state_global->natoms; i++)
840 clear_rvec(state_global->v[i]);
844 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
845 " Ekin, temperature and pressure are incorrect,\n"
846 " the virial will be incorrect when constraints are present.\n"
848 bRerunWarnNoV = FALSE;
852 copy_mat(rerun_fr.box, state_global->box);
853 copy_mat(state_global->box, state->box);
855 if (vsite && (Flags & MD_RERUN_VSITE))
857 if (DOMAINDECOMP(cr))
859 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
863 /* Following is necessary because the graph may get out of sync
864 * with the coordinates if we only have every N'th coordinate set
866 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
867 shift_self(graph, state->box, state->x);
869 construct_vsites(vsite, state->x, ir->delta_t, state->v,
870 top->idef.iparams, top->idef.il,
871 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
874 unshift_self(graph, state->box, state->x);
879 /* Stop Center of Mass motion */
880 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
884 /* for rerun MD always do Neighbour Searching */
885 bNS = (bFirstStep || ir->nstlist != 0);
890 /* Determine whether or not to do Neighbour Searching and LR */
891 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
893 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
894 (ir->nstlist == -1 && nlh.nabnsb > 0));
896 if (bNS && ir->nstlist == -1)
898 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
902 /* check whether we should stop because another simulation has
906 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
907 (multisim_nsteps != ir->nsteps) )
914 "Stopping simulation %d because another one has finished\n",
918 gs.sig[eglsCHKPT] = 1;
923 /* < 0 means stop at next step, > 0 means stop at next NS step */
924 if ( (gs.set[eglsSTOPCOND] < 0) ||
925 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
930 /* Determine whether or not to update the Born radii if doing GB */
931 bBornRadii = bFirstStep;
932 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
937 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
938 do_verbose = bVerbose &&
939 (step % stepout == 0 || bFirstStep || bLastStep);
941 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
949 bMasterState = FALSE;
950 /* Correct the new box if it is too skewed */
951 if (DYNAMIC_BOX(*ir))
953 if (correct_box(fplog, step, state->box, graph))
958 if (DOMAINDECOMP(cr) && bMasterState)
960 dd_collect_state(cr->dd, state, state_global);
964 if (DOMAINDECOMP(cr))
966 /* Repartition the domain decomposition */
967 wallcycle_start(wcycle, ewcDOMDEC);
968 dd_partition_system(fplog, step, cr,
969 bMasterState, nstglobalcomm,
970 state_global, top_global, ir,
971 state, &f, mdatoms, top, fr,
972 vsite, shellfc, constr,
974 do_verbose && !bPMETuneRunning);
975 wallcycle_stop(wcycle, ewcDOMDEC);
976 /* If using an iterative integrator, reallocate space to match the decomposition */
980 if (MASTER(cr) && do_log)
982 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
985 if (ir->efep != efepNO)
987 update_mdatoms(mdatoms, state->lambda[efptMASS]);
990 if ((bRerunMD && rerun_fr.bV) || bExchanged)
993 /* We need the kinetic energy at minus the half step for determining
994 * the full step kinetic energy and possibly for T-coupling.*/
995 /* This may not be quite working correctly yet . . . . */
996 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
997 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
998 constr, NULL, FALSE, state->box,
999 top_global, &bSumEkinhOld,
1000 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1002 clear_mat(force_vir);
1004 /* We write a checkpoint at this MD step when:
1005 * either at an NS step when we signalled through gs,
1006 * or at the last step (but not when we do not want confout),
1007 * but never at the first step or with rerun.
1009 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1010 (bLastStep && (Flags & MD_CONFOUT))) &&
1011 step > ir->init_step && !bRerunMD);
1014 gs.set[eglsCHKPT] = 0;
1017 /* Determine the energy and pressure:
1018 * at nstcalcenergy steps and at energy output steps (set below).
1020 if (EI_VV(ir->eI) && (!bInitStep))
1022 /* for vv, the first half of the integration actually corresponds
1023 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1024 but the virial needs to be calculated on both the current step and the 'next' step. Future
1025 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1027 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1028 bCalcVir = bCalcEner ||
1029 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1033 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1034 bCalcVir = bCalcEner ||
1035 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1038 /* Do we need global communication ? */
1039 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1040 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1041 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1043 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1045 if (do_ene || do_log)
1052 /* these CGLO_ options remain the same throughout the iteration */
1053 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1054 (bGStat ? CGLO_GSTAT : 0)
1057 force_flags = (GMX_FORCE_STATECHANGED |
1058 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1059 GMX_FORCE_ALLFORCES |
1061 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1062 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1063 (bDoFEP ? GMX_FORCE_DHDL : 0)
1068 if (do_per_step(step, ir->nstcalclr))
1070 force_flags |= GMX_FORCE_DO_LR;
1076 /* Now is the time to relax the shells */
1077 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1078 ir, bNS, force_flags,
1081 state, f, force_vir, mdatoms,
1082 nrnb, wcycle, graph, groups,
1083 shellfc, fr, bBornRadii, t, mu_tot,
1095 /* The coordinates (x) are shifted (to get whole molecules)
1097 * This is parallellized as well, and does communication too.
1098 * Check comments in sim_util.c
1100 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1101 state->box, state->x, &state->hist,
1102 f, force_vir, mdatoms, enerd, fcd,
1103 state->lambda, graph,
1104 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1105 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1108 if (bVV && !bStartingFromCpt && !bRerunMD)
1109 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1111 if (ir->eI == eiVV && bInitStep)
1113 /* if using velocity verlet with full time step Ekin,
1114 * take the first half step only to compute the
1115 * virial for the first step. From there,
1116 * revert back to the initial coordinates
1117 * so that the input is actually the initial step.
1119 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1123 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1124 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1127 /* If we are using twin-range interactions where the long-range component
1128 * is only evaluated every nstcalclr>1 steps, we should do a special update
1129 * step to combine the long-range forces on these steps.
1130 * For nstcalclr=1 this is not done, since the forces would have been added
1131 * directly to the short-range forces already.
1133 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1135 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1136 f, bUpdateDoLR, fr->f_twin, fcd,
1137 ekind, M, upd, bInitStep, etrtVELOCITY1,
1138 cr, nrnb, constr, &top->idef);
1140 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1142 gmx_iterate_init(&iterate, TRUE);
1144 /* for iterations, we save these vectors, as we will be self-consistently iterating
1147 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1149 /* save the state */
1150 if (iterate.bIterationActive)
1152 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1155 bFirstIterate = TRUE;
1156 while (bFirstIterate || iterate.bIterationActive)
1158 if (iterate.bIterationActive)
1160 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1161 if (bFirstIterate && bTrotter)
1163 /* The first time through, we need a decent first estimate
1164 of veta(t+dt) to compute the constraints. Do
1165 this by computing the box volume part of the
1166 trotter integration at this time. Nothing else
1167 should be changed by this routine here. If
1168 !(first time), we start with the previous value
1171 veta_save = state->veta;
1172 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1173 vetanew = state->veta;
1174 state->veta = veta_save;
1179 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1181 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1182 state, fr->bMolPBC, graph, f,
1183 &top->idef, shake_vir,
1184 cr, nrnb, wcycle, upd, constr,
1185 TRUE, bCalcVir, vetanew);
1189 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1195 /* Need to unshift here if a do_force has been
1196 called in the previous step */
1197 unshift_self(graph, state->box, state->x);
1200 /* if VV, compute the pressure and constraints */
1201 /* For VV2, we strictly only need this if using pressure
1202 * control, but we really would like to have accurate pressures
1204 * Think about ways around this in the future?
1205 * For now, keep this choice in comments.
1207 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1208 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1210 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1211 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1213 bSumEkinhOld = TRUE;
1215 /* for vv, the first half of the integration actually corresponds to the previous step.
1216 So we need information from the last step in the first half of the integration */
1217 if (bGStat || do_per_step(step-1, nstglobalcomm))
1219 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1220 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1221 constr, NULL, FALSE, state->box,
1222 top_global, &bSumEkinhOld,
1225 | (bTemp ? CGLO_TEMPERATURE : 0)
1226 | (bPres ? CGLO_PRESSURE : 0)
1227 | (bPres ? CGLO_CONSTRAINT : 0)
1228 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1229 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1232 /* explanation of above:
1233 a) We compute Ekin at the full time step
1234 if 1) we are using the AveVel Ekin, and it's not the
1235 initial step, or 2) if we are using AveEkin, but need the full
1236 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1237 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1238 EkinAveVel because it's needed for the pressure */
1240 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1245 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1246 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1253 /* We need the kinetic energy at minus the half step for determining
1254 * the full step kinetic energy and possibly for T-coupling.*/
1255 /* This may not be quite working correctly yet . . . . */
1256 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1257 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1258 constr, NULL, FALSE, state->box,
1259 top_global, &bSumEkinhOld,
1260 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1265 if (iterate.bIterationActive &&
1266 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1267 state->veta, &vetanew))
1271 bFirstIterate = FALSE;
1274 if (bTrotter && !bInitStep)
1276 copy_mat(shake_vir, state->svir_prev);
1277 copy_mat(force_vir, state->fvir_prev);
1278 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1280 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1281 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1282 enerd->term[F_EKIN] = trace(ekind->ekin);
1285 /* if it's the initial step, we performed this first step just to get the constraint virial */
1286 if (bInitStep && ir->eI == eiVV)
1288 copy_rvecn(cbuf, state->v, 0, state->natoms);
1292 /* MRS -- now done iterating -- compute the conserved quantity */
1295 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1298 last_ekin = enerd->term[F_EKIN];
1300 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1302 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1304 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1307 sum_dhdl(enerd, state->lambda, ir->fepvals);
1311 /* ######## END FIRST UPDATE STEP ############## */
1312 /* ######## If doing VV, we now have v(dt) ###### */
1315 /* perform extended ensemble sampling in lambda - we don't
1316 actually move to the new state before outputting
1317 statistics, but if performing simulated tempering, we
1318 do update the velocities and the tau_t. */
1320 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, mcrng, state->v, mdatoms);
1321 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1322 copy_df_history(&state_global->dfhist,&state->dfhist);
1325 /* Now we have the energies and forces corresponding to the
1326 * coordinates at time t. We must output all of this before
1329 do_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1330 ir, state, state_global, top_global, fr, upd,
1331 outf, mdebin, ekind, f, f_global,
1332 wcycle, mcrng, &nchkpt,
1333 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1336 /* kludge -- virial is lost with restart for NPT control. Must restart */
1337 if (bStartingFromCpt && bVV)
1339 copy_mat(state->svir_prev, shake_vir);
1340 copy_mat(state->fvir_prev, force_vir);
1343 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1345 /* Check whether everything is still allright */
1346 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1347 #ifdef GMX_THREAD_MPI
1352 /* this is just make gs.sig compatible with the hack
1353 of sending signals around by MPI_Reduce with together with
1355 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1357 gs.sig[eglsSTOPCOND] = 1;
1359 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1361 gs.sig[eglsSTOPCOND] = -1;
1363 /* < 0 means stop at next step, > 0 means stop at next NS step */
1367 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1368 gmx_get_signal_name(),
1369 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1373 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1374 gmx_get_signal_name(),
1375 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1377 handled_stop_condition = (int)gmx_get_stop_condition();
1379 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1380 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1381 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1383 /* Signal to terminate the run */
1384 gs.sig[eglsSTOPCOND] = 1;
1387 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1389 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1392 if (bResetCountersHalfMaxH && MASTER(cr) &&
1393 elapsed_time > max_hours*60.0*60.0*0.495)
1395 gs.sig[eglsRESETCOUNTERS] = 1;
1398 if (ir->nstlist == -1 && !bRerunMD)
1400 /* When bGStatEveryStep=FALSE, global_stat is only called
1401 * when we check the atom displacements, not at NS steps.
1402 * This means that also the bonded interaction count check is not
1403 * performed immediately after NS. Therefore a few MD steps could
1404 * be performed with missing interactions.
1405 * But wrong energies are never written to file,
1406 * since energies are only written after global_stat
1409 if (step >= nlh.step_nscheck)
1411 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1412 nlh.scale_tot, state->x);
1416 /* This is not necessarily true,
1417 * but step_nscheck is determined quite conservatively.
1423 /* In parallel we only have to check for checkpointing in steps
1424 * where we do global communication,
1425 * otherwise the other nodes don't know.
1427 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1430 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1431 gs.set[eglsCHKPT] == 0)
1433 gs.sig[eglsCHKPT] = 1;
1436 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1441 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1443 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1445 gmx_bool bIfRandomize;
1446 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1447 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1448 if (constr && bIfRandomize)
1450 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1451 state, fr->bMolPBC, graph, f,
1452 &top->idef, tmp_vir,
1453 cr, nrnb, wcycle, upd, constr,
1454 TRUE, bCalcVir, vetanew);
1459 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1461 gmx_iterate_init(&iterate, TRUE);
1462 /* for iterations, we save these vectors, as we will be redoing the calculations */
1463 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1466 bFirstIterate = TRUE;
1467 while (bFirstIterate || iterate.bIterationActive)
1469 /* We now restore these vectors to redo the calculation with improved extended variables */
1470 if (iterate.bIterationActive)
1472 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1475 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1476 so scroll down for that logic */
1478 /* ######### START SECOND UPDATE STEP ################# */
1479 /* Box is changed in update() when we do pressure coupling,
1480 * but we should still use the old box for energy corrections and when
1481 * writing it to the energy file, so it matches the trajectory files for
1482 * the same timestep above. Make a copy in a separate array.
1484 copy_mat(state->box, lastbox);
1489 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1491 wallcycle_start(wcycle, ewcUPDATE);
1492 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1495 if (iterate.bIterationActive)
1503 /* we use a new value of scalevir to converge the iterations faster */
1504 scalevir = tracevir/trace(shake_vir);
1506 msmul(shake_vir, scalevir, shake_vir);
1507 m_add(force_vir, shake_vir, total_vir);
1508 clear_mat(shake_vir);
1510 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1511 /* We can only do Berendsen coupling after we have summed
1512 * the kinetic energy or virial. Since the happens
1513 * in global_state after update, we should only do it at
1514 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1519 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1520 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1525 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1527 /* velocity half-step update */
1528 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1529 bUpdateDoLR, fr->f_twin, fcd,
1530 ekind, M, upd, FALSE, etrtVELOCITY2,
1531 cr, nrnb, constr, &top->idef);
1534 /* Above, initialize just copies ekinh into ekin,
1535 * it doesn't copy position (for VV),
1536 * and entire integrator for MD.
1539 if (ir->eI == eiVVAK)
1541 copy_rvecn(state->x, cbuf, 0, state->natoms);
1543 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1545 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1546 bUpdateDoLR, fr->f_twin, fcd,
1547 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1548 wallcycle_stop(wcycle, ewcUPDATE);
1550 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1551 fr->bMolPBC, graph, f,
1552 &top->idef, shake_vir,
1553 cr, nrnb, wcycle, upd, constr,
1554 FALSE, bCalcVir, state->veta);
1556 if (ir->eI == eiVVAK)
1558 /* erase F_EKIN and F_TEMP here? */
1559 /* just compute the kinetic energy at the half step to perform a trotter step */
1560 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1561 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1562 constr, NULL, FALSE, lastbox,
1563 top_global, &bSumEkinhOld,
1564 cglo_flags | CGLO_TEMPERATURE
1566 wallcycle_start(wcycle, ewcUPDATE);
1567 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1568 /* now we know the scaling, we can compute the positions again again */
1569 copy_rvecn(cbuf, state->x, 0, state->natoms);
1571 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1573 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1574 bUpdateDoLR, fr->f_twin, fcd,
1575 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1576 wallcycle_stop(wcycle, ewcUPDATE);
1578 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1579 /* are the small terms in the shake_vir here due
1580 * to numerical errors, or are they important
1581 * physically? I'm thinking they are just errors, but not completely sure.
1582 * For now, will call without actually constraining, constr=NULL*/
1583 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1584 state, fr->bMolPBC, graph, f,
1585 &top->idef, tmp_vir,
1586 cr, nrnb, wcycle, upd, NULL,
1592 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1595 if (fr->bSepDVDL && fplog && do_log)
1597 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1601 /* this factor or 2 correction is necessary
1602 because half of the constraint force is removed
1603 in the vv step, so we have to double it. See
1604 the Redmine issue #1255. It is not yet clear
1605 if the factor of 2 is exact, or just a very
1606 good approximation, and this will be
1607 investigated. The next step is to see if this
1608 can be done adding a dhdl contribution from the
1609 rattle step, but this is somewhat more
1610 complicated with the current code. Will be
1611 investigated, hopefully for 4.6.3. However,
1612 this current solution is much better than
1613 having it completely wrong.
1615 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1619 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1624 /* Need to unshift here */
1625 unshift_self(graph, state->box, state->x);
1630 wallcycle_start(wcycle, ewcVSITECONSTR);
1633 shift_self(graph, state->box, state->x);
1635 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1636 top->idef.iparams, top->idef.il,
1637 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1641 unshift_self(graph, state->box, state->x);
1643 wallcycle_stop(wcycle, ewcVSITECONSTR);
1646 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1647 /* With Leap-Frog we can skip compute_globals at
1648 * non-communication steps, but we need to calculate
1649 * the kinetic energy one step before communication.
1651 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1653 if (ir->nstlist == -1 && bFirstIterate)
1655 gs.sig[eglsNABNSB] = nlh.nabnsb;
1657 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1658 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1660 bFirstIterate ? &gs : NULL,
1661 (step_rel % gs.nstms == 0) &&
1662 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1664 top_global, &bSumEkinhOld,
1666 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1667 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1668 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1669 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1670 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1671 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1674 if (ir->nstlist == -1 && bFirstIterate)
1676 nlh.nabnsb = gs.set[eglsNABNSB];
1677 gs.set[eglsNABNSB] = 0;
1680 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1681 /* ############# END CALC EKIN AND PRESSURE ################# */
1683 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1684 the virial that should probably be addressed eventually. state->veta has better properies,
1685 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1686 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1688 if (iterate.bIterationActive &&
1689 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1690 trace(shake_vir), &tracevir))
1694 bFirstIterate = FALSE;
1697 if (!bVV || bRerunMD)
1699 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1700 sum_dhdl(enerd, state->lambda, ir->fepvals);
1702 update_box(fplog, step, ir, mdatoms, state, f,
1703 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1705 /* ################# END UPDATE STEP 2 ################# */
1706 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1708 /* The coordinates (x) were unshifted in update */
1711 /* We will not sum ekinh_old,
1712 * so signal that we still have to do it.
1714 bSumEkinhOld = TRUE;
1717 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1719 /* use the directly determined last velocity, not actually the averaged half steps */
1720 if (bTrotter && ir->eI == eiVV)
1722 enerd->term[F_EKIN] = last_ekin;
1724 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1728 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1732 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1734 /* ######### END PREPARING EDR OUTPUT ########### */
1739 gmx_bool do_dr, do_or;
1741 if (fplog && do_log && bDoExpanded)
1743 /* only needed if doing expanded ensemble */
1744 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1745 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1747 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1751 upd_mdebin(mdebin, bDoDHDL, TRUE,
1752 t, mdatoms->tmass, enerd, state,
1753 ir->fepvals, ir->expandedvals, lastbox,
1754 shake_vir, force_vir, total_vir, pres,
1755 ekind, mu_tot, constr);
1759 upd_mdebin_step(mdebin);
1762 do_dr = do_per_step(step, ir->nstdisreout);
1763 do_or = do_per_step(step, ir->nstorireout);
1765 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
1767 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1769 if (ir->ePull != epullNO)
1771 pull_print_output(ir->pull, step, t);
1774 if (do_per_step(step, ir->nstlog))
1776 if (fflush(fplog) != 0)
1778 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1784 /* Have to do this part _after_ outputting the logfile and the edr file */
1785 /* Gets written into the state at the beginning of next loop*/
1786 state->fep_state = lamnew;
1788 /* Print the remaining wall clock time for the run */
1789 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1793 fprintf(stderr, "\n");
1795 print_time(stderr, walltime_accounting, step, ir, cr);
1798 /* Replica exchange */
1800 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1801 do_per_step(step, repl_ex_nst))
1803 bExchanged = replica_exchange(fplog, cr, repl_ex,
1804 state_global, enerd,
1807 if (bExchanged && DOMAINDECOMP(cr))
1809 dd_partition_system(fplog, step, cr, TRUE, 1,
1810 state_global, top_global, ir,
1811 state, &f, mdatoms, top, fr,
1812 vsite, shellfc, constr,
1813 nrnb, wcycle, FALSE);
1819 bStartingFromCpt = FALSE;
1821 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1822 /* With all integrators, except VV, we need to retain the pressure
1823 * at the current step for coupling at the next step.
1825 if ((state->flags & (1<<estPRES_PREV)) &&
1827 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1829 /* Store the pressure in t_state for pressure coupling
1830 * at the next MD step.
1832 copy_mat(pres, state->pres_prev);
1835 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1837 if ( (membed != NULL) && (!bLastStep) )
1839 rescale_membed(step_rel, membed, state_global->x);
1846 /* read next frame from input trajectory */
1847 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1852 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1856 if (!bRerunMD || !rerun_fr.bStep)
1858 /* increase the MD step number */
1863 cycles = wallcycle_stop(wcycle, ewcSTEP);
1864 if (DOMAINDECOMP(cr) && wcycle)
1866 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1869 if (bPMETuneRunning || bPMETuneTry)
1871 /* PME grid + cut-off optimization with GPUs or PME nodes */
1873 /* Count the total cycles over the last steps */
1874 cycles_pmes += cycles;
1876 /* We can only switch cut-off at NS steps */
1877 if (step % ir->nstlist == 0)
1879 /* PME grid + cut-off optimization with GPUs or PME nodes */
1882 if (DDMASTER(cr->dd))
1884 /* PME node load is too high, start tuning */
1885 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1887 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1889 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1891 bPMETuneTry = FALSE;
1894 if (bPMETuneRunning)
1896 /* init_step might not be a multiple of nstlist,
1897 * but the first cycle is always skipped anyhow.
1900 pme_load_balance(pme_loadbal, cr,
1901 (bVerbose && MASTER(cr)) ? stderr : NULL,
1903 ir, state, cycles_pmes,
1904 fr->ic, fr->nbv, &fr->pmedata,
1907 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1908 fr->ewaldcoeff = fr->ic->ewaldcoeff;
1909 fr->rlist = fr->ic->rlist;
1910 fr->rlistlong = fr->ic->rlistlong;
1911 fr->rcoulomb = fr->ic->rcoulomb;
1912 fr->rvdw = fr->ic->rvdw;
1918 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1919 gs.set[eglsRESETCOUNTERS] != 0)
1921 /* Reset all the counters related to performance over the run */
1922 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1923 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1924 wcycle_set_reset_counters(wcycle, -1);
1925 if (!(cr->duty & DUTY_PME))
1927 /* Tell our PME node to reset its counters */
1928 gmx_pme_send_resetcounters(cr, step);
1930 /* Correct max_hours for the elapsed time */
1931 max_hours -= elapsed_time/(60.0*60.0);
1932 bResetCountersHalfMaxH = FALSE;
1933 gs.set[eglsRESETCOUNTERS] = 0;
1937 /* End of main MD loop */
1940 /* Stop measuring walltime */
1941 walltime_accounting_end(walltime_accounting);
1943 if (bRerunMD && MASTER(cr))
1948 if (!(cr->duty & DUTY_PME))
1950 /* Tell the PME only node to finish */
1951 gmx_pme_send_finish(cr);
1956 if (ir->nstcalcenergy > 0 && !bRerunMD)
1958 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
1959 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1967 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1969 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)));
1970 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1973 if (pme_loadbal != NULL)
1975 pme_loadbal_done(pme_loadbal, cr, fplog,
1976 fr->nbv != NULL && fr->nbv->bUseGPU);
1979 if (shellfc && fplog)
1981 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1982 (nconverged*100.0)/step_rel);
1983 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1987 if (repl_ex_nst > 0 && MASTER(cr))
1989 print_replica_exchange_statistics(fplog, repl_ex);
1992 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);