1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
56 #include "md_support.h"
57 #include "md_logging.h"
73 #include "domdec_network.h"
79 #include "compute_io.h"
81 #include "checkpoint.h"
82 #include "mtop_util.h"
83 #include "sighandler.h"
86 #include "pme_loadbal.h"
89 #include "types/nlistheuristics.h"
90 #include "types/iteratedconstraints.h"
91 #include "nbnxn_cuda_data_mgmt.h"
93 #include "gromacs/utility/gmxmpi.h"
99 static void reset_all_counters(FILE *fplog, t_commrec *cr,
100 gmx_large_int_t step,
101 gmx_large_int_t *step_rel, t_inputrec *ir,
102 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
103 gmx_runtime_t *runtime,
104 nbnxn_cuda_ptr_t cu_nbv)
106 char sbuf[STEPSTRSIZE];
108 /* Reset all the counters related to performance over the run */
109 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
110 gmx_step_str(step, sbuf));
114 nbnxn_cuda_reset_timings(cu_nbv);
117 wallcycle_stop(wcycle, ewcRUN);
118 wallcycle_reset_all(wcycle);
119 if (DOMAINDECOMP(cr))
121 reset_dd_statistics_counters(cr->dd);
124 ir->init_step += *step_rel;
125 ir->nsteps -= *step_rel;
127 wallcycle_start(wcycle, ewcRUN);
128 runtime_start(runtime);
129 print_date_and_time(fplog, cr->nodeid, "Restarted time", runtime);
132 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
133 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
135 gmx_vsite_t *vsite, gmx_constr_t constr,
136 int stepout, t_inputrec *ir,
137 gmx_mtop_t *top_global,
139 t_state *state_global,
141 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
142 gmx_edsam_t ed, t_forcerec *fr,
143 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
144 real cpt_period, real max_hours,
145 const char gmx_unused *deviceOptions,
147 gmx_runtime_t *runtime)
150 gmx_large_int_t step, step_rel;
152 double t, t0, lam0[efptNR];
153 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
154 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
155 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
156 bBornRadii, bStartingFromCpt;
157 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
158 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
159 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 df_history_t df_history;
177 t_state *state = NULL;
178 rvec *f_global = NULL;
181 gmx_enerdata_t *enerd;
183 gmx_global_stat_t gstat;
184 gmx_update_t upd = NULL;
185 t_graph *graph = NULL;
187 gmx_rng_t mcrng = NULL;
189 gmx_groups_t *groups;
190 gmx_ekindata_t *ekind, *ekind_save;
191 gmx_shellfc_t shellfc;
192 int count, nconverged = 0;
195 gmx_bool bIonize = FALSE;
196 gmx_bool bTCR = FALSE, bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
198 gmx_bool bResetCountersHalfMaxH = FALSE;
199 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
200 gmx_bool bUpdateDoLR;
201 real mu_aver = 0, dvdl_constr;
202 int a0, a1, gnx = 0, ii;
203 atom_id *grpindex = NULL;
205 t_coupl_rec *tcr = NULL;
206 rvec *xcopy = NULL, *vcopy = NULL, *cbuf = NULL;
207 matrix boxcopy = {{0}}, lastbox;
209 real fom, oldfom, veta_save, pcurr, scalevir, tracevir;
216 real saved_conserved_quantity = 0;
221 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
222 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
223 gmx_iterate_t iterate;
224 gmx_large_int_t multisim_nsteps = -1; /* number of steps to do before first multisim
225 simulation stops. If equal to zero, don't
226 communicate any more between multisims.*/
227 /* PME load balancing data for GPU kernels */
228 pme_load_balancing_t pme_loadbal = NULL;
230 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
233 /* Temporary addition for FAHCORE checkpointing */
237 /* Check for special mdrun options */
238 bRerunMD = (Flags & MD_RERUN);
239 bIonize = (Flags & MD_IONIZE);
240 bFFscan = (Flags & MD_FFSCAN);
241 bAppend = (Flags & MD_APPENDFILES);
242 if (Flags & MD_RESETCOUNTERSHALFWAY)
246 /* Signal to reset the counters half the simulation steps. */
247 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
249 /* Signal to reset the counters halfway the simulation time. */
250 bResetCountersHalfMaxH = (max_hours > 0);
253 /* md-vv uses averaged full step velocities for T-control
254 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
255 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
257 if (bVV) /* to store the initial velocities while computing virial */
259 snew(cbuf, top_global->natoms);
261 /* all the iteratative cases - only if there are constraints */
262 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
263 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
264 false in this step. The correct value, true or false,
265 is set at each step, as it depends on the frequency of temperature
266 and pressure control.*/
267 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
271 /* Since we don't know if the frames read are related in any way,
272 * rebuild the neighborlist at every step.
275 ir->nstcalcenergy = 1;
279 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
281 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
282 bGStatEveryStep = (nstglobalcomm == 1);
284 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
287 "To reduce the energy communication with nstlist = -1\n"
288 "the neighbor list validity should not be checked at every step,\n"
289 "this means that exact integration is not guaranteed.\n"
290 "The neighbor list validity is checked after:\n"
291 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
292 "In most cases this will result in exact integration.\n"
293 "This reduces the energy communication by a factor of 2 to 3.\n"
294 "If you want less energy communication, set nstlist > 3.\n\n");
297 if (bRerunMD || bFFscan)
301 groups = &top_global->groups;
304 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
305 &(state_global->fep_state), lam0,
306 nrnb, top_global, &upd,
307 nfile, fnm, &outf, &mdebin,
308 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
310 clear_mat(total_vir);
312 /* Energy terms and groups */
314 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
316 if (DOMAINDECOMP(cr))
322 snew(f, top_global->natoms);
325 /* lambda Monte carlo random number generator */
328 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
330 /* copy the state into df_history */
331 copy_df_history(&df_history, &state_global->dfhist);
333 /* Kinetic energy data */
335 init_ekindata(fplog, top_global, &(ir->opts), ekind);
336 /* needed for iteration of constraints */
338 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
339 /* Copy the cos acceleration to the groups struct */
340 ekind->cosacc.cos_accel = ir->cos_accel;
342 gstat = global_stat_init(ir);
345 /* Check for polarizable models and flexible constraints */
346 shellfc = init_shell_flexcon(fplog,
347 top_global, n_flexible_constraints(constr),
348 (ir->bContinuation ||
349 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
350 NULL : state_global->x);
354 #ifdef GMX_THREAD_MPI
355 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
357 set_deform_reference_box(upd,
358 deform_init_init_step_tpx,
359 deform_init_box_tpx);
360 #ifdef GMX_THREAD_MPI
361 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
366 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
367 if ((io > 2000) && MASTER(cr))
370 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
375 if (DOMAINDECOMP(cr))
377 top = dd_init_local_top(top_global);
380 dd_init_local_state(cr->dd, state_global, state);
382 if (DDMASTER(cr->dd) && ir->nstfout)
384 snew(f_global, state_global->natoms);
391 /* Initialize the particle decomposition and split the topology */
392 top = split_system(fplog, top_global, ir, cr);
394 pd_cg_range(cr, &fr->cg0, &fr->hcg);
395 pd_at_range(cr, &a0, &a1);
399 top = gmx_mtop_generate_local_top(top_global, ir);
402 a1 = top_global->natoms;
405 forcerec_set_excl_load(fr, top, cr);
407 state = partdec_init_local_state(cr, state_global);
410 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
414 set_vsite_top(vsite, top, mdatoms, cr);
417 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
419 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
424 make_local_shells(cr, mdatoms, shellfc);
427 init_bonded_thread_force_reduction(fr, &top->idef);
429 if (ir->pull && PAR(cr))
431 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
435 if (DOMAINDECOMP(cr))
437 /* Distribute the charge groups over the nodes from the master node */
438 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
439 state_global, top_global, ir,
440 state, &f, mdatoms, top, fr,
441 vsite, shellfc, constr,
442 nrnb, wcycle, FALSE);
446 update_mdatoms(mdatoms, state->lambda[efptMASS]);
448 if (opt2bSet("-cpi", nfile, fnm))
450 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
454 bStateFromCP = FALSE;
461 /* Update mdebin with energy history if appending to output files */
462 if (Flags & MD_APPENDFILES)
464 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
468 /* We might have read an energy history from checkpoint,
469 * free the allocated memory and reset the counts.
471 done_energyhistory(&state_global->enerhist);
472 init_energyhistory(&state_global->enerhist);
475 /* Set the initial energy history in state by updating once */
476 update_energyhistory(&state_global->enerhist, mdebin);
479 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
481 /* Set the random state if we read a checkpoint file */
482 set_stochd_state(upd, state);
485 if (state->flags & (1<<estMC_RNG))
487 set_mc_state(mcrng, state);
490 /* Initialize constraints */
493 if (!DOMAINDECOMP(cr))
495 set_constraints(constr, top, ir, mdatoms, cr);
499 /* Check whether we have to GCT stuff */
500 bTCR = ftp2bSet(efGCT, nfile, fnm);
505 fprintf(stderr, "Will do General Coupling Theory!\n");
507 gnx = top_global->mols.nr;
509 for (i = 0; (i < gnx); i++)
517 /* We need to be sure replica exchange can only occur
518 * when the energies are current */
519 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
520 "repl_ex_nst", &repl_ex_nst);
521 /* This check needs to happen before inter-simulation
522 * signals are initialized, too */
524 if (repl_ex_nst > 0 && MASTER(cr))
526 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
527 repl_ex_nst, repl_ex_nex, repl_ex_seed);
530 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
531 * With perturbed charges with soft-core we should not change the cut-off.
533 if ((Flags & MD_TUNEPME) &&
534 EEL_PME(fr->eeltype) &&
535 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
536 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
539 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
541 if (cr->duty & DUTY_PME)
543 /* Start tuning right away, as we can't measure the load */
544 bPMETuneRunning = TRUE;
548 /* Separate PME nodes, we can measure the PP/PME load balance */
553 if (!ir->bContinuation && !bRerunMD)
555 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
557 /* Set the velocities of frozen particles to zero */
558 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
560 for (m = 0; m < DIM; m++)
562 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
572 /* Constrain the initial coordinates and velocities */
573 do_constrain_first(fplog, constr, ir, mdatoms, state,
578 /* Construct the virtual sites for the initial configuration */
579 construct_vsites(vsite, state->x, ir->delta_t, NULL,
580 top->idef.iparams, top->idef.il,
581 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
587 /* set free energy calculation frequency as the minimum
588 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
589 nstfep = ir->fepvals->nstdhdl;
592 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
596 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
599 /* I'm assuming we need global communication the first time! MRS */
600 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
601 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
602 | (bVV ? CGLO_PRESSURE : 0)
603 | (bVV ? CGLO_CONSTRAINT : 0)
604 | (bRerunMD ? CGLO_RERUNMD : 0)
605 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
607 bSumEkinhOld = FALSE;
608 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
609 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
610 constr, NULL, FALSE, state->box,
611 top_global, &pcurr, &bSumEkinhOld, cglo_flags);
612 if (ir->eI == eiVVAK)
614 /* a second call to get the half step temperature initialized as well */
615 /* we do the same call as above, but turn the pressure off -- internally to
616 compute_globals, this is recognized as a velocity verlet half-step
617 kinetic energy calculation. This minimized excess variables, but
618 perhaps loses some logic?*/
620 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
621 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
622 constr, NULL, FALSE, state->box,
623 top_global, &pcurr, &bSumEkinhOld,
624 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
627 /* Calculate the initial half step temperature, and save the ekinh_old */
628 if (!(Flags & MD_STARTFROMCPT))
630 for (i = 0; (i < ir->opts.ngtc); i++)
632 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
637 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
638 and there is no previous step */
641 /* if using an iterative algorithm, we need to create a working directory for the state. */
644 bufstate = init_bufstate(state);
648 snew(xcopy, state->natoms);
649 snew(vcopy, state->natoms);
650 copy_rvecn(state->x, xcopy, 0, state->natoms);
651 copy_rvecn(state->v, vcopy, 0, state->natoms);
652 copy_mat(state->box, boxcopy);
655 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
656 temperature control */
657 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
661 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
664 "RMS relative constraint deviation after constraining: %.2e\n",
665 constr_rmsd(constr, FALSE));
667 if (EI_STATE_VELOCITY(ir->eI))
669 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
673 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
674 " input trajectory '%s'\n\n",
675 *(top_global->name), opt2fn("-rerun", nfile, fnm));
678 fprintf(stderr, "Calculated time to finish depends on nsteps from "
679 "run input file,\nwhich may not correspond to the time "
680 "needed to process input trajectory.\n\n");
686 fprintf(stderr, "starting mdrun '%s'\n",
687 *(top_global->name));
690 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
694 sprintf(tbuf, "%s", "infinite");
696 if (ir->init_step > 0)
698 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
699 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
700 gmx_step_str(ir->init_step, sbuf2),
701 ir->init_step*ir->delta_t);
705 fprintf(stderr, "%s steps, %s ps.\n",
706 gmx_step_str(ir->nsteps, sbuf), tbuf);
709 fprintf(fplog, "\n");
712 /* Set and write start time */
713 runtime_start(runtime);
714 print_date_and_time(fplog, cr->nodeid, "Started mdrun", runtime);
715 wallcycle_start(wcycle, ewcRUN);
718 fprintf(fplog, "\n");
721 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
723 chkpt_ret = fcCheckPointParallel( cr->nodeid,
727 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
732 /***********************************************************
736 ************************************************************/
738 /* if rerunMD then read coordinates and velocities from input trajectory */
741 if (getenv("GMX_FORCE_UPDATE"))
749 bNotLastFrame = read_first_frame(oenv, &status,
750 opt2fn("-rerun", nfile, fnm),
751 &rerun_fr, TRX_NEED_X | TRX_READ_V);
752 if (rerun_fr.natoms != top_global->natoms)
755 "Number of atoms in trajectory (%d) does not match the "
756 "run input file (%d)\n",
757 rerun_fr.natoms, top_global->natoms);
759 if (ir->ePBC != epbcNONE)
763 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);
765 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
767 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
774 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
777 if (ir->ePBC != epbcNONE)
779 /* Set the shift vectors.
780 * Necessary here when have a static box different from the tpr box.
782 calc_shifts(rerun_fr.box, fr->shift_vec);
786 /* loop over MD steps or if rerunMD to end of input trajectory */
788 /* Skip the first Nose-Hoover integration when we get the state from tpx */
789 bStateFromTPX = !bStateFromCP;
790 bInitStep = bFirstStep && (bStateFromTPX || bVV);
791 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
793 bSumEkinhOld = FALSE;
796 init_global_signals(&gs, cr, ir, repl_ex_nst);
798 step = ir->init_step;
801 if (ir->nstlist == -1)
803 init_nlistheuristics(&nlh, bGStatEveryStep, step);
806 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
808 /* check how many steps are left in other sims */
809 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
813 /* and stop now if we should */
814 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
815 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
816 while (!bLastStep || (bRerunMD && bNotLastFrame))
819 wallcycle_start(wcycle, ewcSTEP);
825 step = rerun_fr.step;
826 step_rel = step - ir->init_step;
839 bLastStep = (step_rel == ir->nsteps);
840 t = t0 + step*ir->delta_t;
843 if (ir->efep != efepNO || ir->bSimTemp)
845 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
846 requiring different logic. */
848 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
849 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
850 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
851 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
856 update_annealing_target_temp(&(ir->opts), t);
861 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
863 for (i = 0; i < state_global->natoms; i++)
865 copy_rvec(rerun_fr.x[i], state_global->x[i]);
869 for (i = 0; i < state_global->natoms; i++)
871 copy_rvec(rerun_fr.v[i], state_global->v[i]);
876 for (i = 0; i < state_global->natoms; i++)
878 clear_rvec(state_global->v[i]);
882 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
883 " Ekin, temperature and pressure are incorrect,\n"
884 " the virial will be incorrect when constraints are present.\n"
886 bRerunWarnNoV = FALSE;
890 copy_mat(rerun_fr.box, state_global->box);
891 copy_mat(state_global->box, state->box);
893 if (vsite && (Flags & MD_RERUN_VSITE))
895 if (DOMAINDECOMP(cr))
897 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
901 /* Following is necessary because the graph may get out of sync
902 * with the coordinates if we only have every N'th coordinate set
904 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
905 shift_self(graph, state->box, state->x);
907 construct_vsites(vsite, state->x, ir->delta_t, state->v,
908 top->idef.iparams, top->idef.il,
909 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
912 unshift_self(graph, state->box, state->x);
917 /* Stop Center of Mass motion */
918 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
920 /* Copy back starting coordinates in case we're doing a forcefield scan */
923 for (ii = 0; (ii < state->natoms); ii++)
925 copy_rvec(xcopy[ii], state->x[ii]);
926 copy_rvec(vcopy[ii], state->v[ii]);
928 copy_mat(boxcopy, state->box);
933 /* for rerun MD always do Neighbour Searching */
934 bNS = (bFirstStep || ir->nstlist != 0);
939 /* Determine whether or not to do Neighbour Searching and LR */
940 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
942 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
943 (ir->nstlist == -1 && nlh.nabnsb > 0));
945 if (bNS && ir->nstlist == -1)
947 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
951 /* check whether we should stop because another simulation has
955 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
956 (multisim_nsteps != ir->nsteps) )
963 "Stopping simulation %d because another one has finished\n",
967 gs.sig[eglsCHKPT] = 1;
972 /* < 0 means stop at next step, > 0 means stop at next NS step */
973 if ( (gs.set[eglsSTOPCOND] < 0) ||
974 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
979 /* Determine whether or not to update the Born radii if doing GB */
980 bBornRadii = bFirstStep;
981 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
986 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
987 do_verbose = bVerbose &&
988 (step % stepout == 0 || bFirstStep || bLastStep);
990 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
998 bMasterState = FALSE;
999 /* Correct the new box if it is too skewed */
1000 if (DYNAMIC_BOX(*ir))
1002 if (correct_box(fplog, step, state->box, graph))
1004 bMasterState = TRUE;
1007 if (DOMAINDECOMP(cr) && bMasterState)
1009 dd_collect_state(cr->dd, state, state_global);
1013 if (DOMAINDECOMP(cr))
1015 /* Repartition the domain decomposition */
1016 wallcycle_start(wcycle, ewcDOMDEC);
1017 dd_partition_system(fplog, step, cr,
1018 bMasterState, nstglobalcomm,
1019 state_global, top_global, ir,
1020 state, &f, mdatoms, top, fr,
1021 vsite, shellfc, constr,
1023 do_verbose && !bPMETuneRunning);
1024 wallcycle_stop(wcycle, ewcDOMDEC);
1025 /* If using an iterative integrator, reallocate space to match the decomposition */
1029 if (MASTER(cr) && do_log && !bFFscan)
1031 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
1034 if (ir->efep != efepNO)
1036 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1039 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1042 /* We need the kinetic energy at minus the half step for determining
1043 * the full step kinetic energy and possibly for T-coupling.*/
1044 /* This may not be quite working correctly yet . . . . */
1045 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1046 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1047 constr, NULL, FALSE, state->box,
1048 top_global, &pcurr, &bSumEkinhOld,
1049 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1051 clear_mat(force_vir);
1053 /* Ionize the atoms if necessary */
1056 ionize(fplog, oenv, mdatoms, top_global, t, ir, state->x, state->v,
1057 mdatoms->start, mdatoms->start+mdatoms->homenr, state->box, cr);
1060 /* Update force field in ffscan program */
1063 if (update_forcefield(fplog,
1065 mdatoms->nr, state->x, state->box))
1073 /* We write a checkpoint at this MD step when:
1074 * either at an NS step when we signalled through gs,
1075 * or at the last step (but not when we do not want confout),
1076 * but never at the first step or with rerun.
1078 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1079 (bLastStep && (Flags & MD_CONFOUT))) &&
1080 step > ir->init_step && !bRerunMD);
1083 gs.set[eglsCHKPT] = 0;
1086 /* Determine the energy and pressure:
1087 * at nstcalcenergy steps and at energy output steps (set below).
1089 if (EI_VV(ir->eI) && (!bInitStep))
1091 /* for vv, the first half of the integration actually corresponds
1092 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1093 but the virial needs to be calculated on both the current step and the 'next' step. Future
1094 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1096 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1097 bCalcVir = bCalcEner ||
1098 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1102 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1103 bCalcVir = bCalcEner ||
1104 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1107 /* Do we need global communication ? */
1108 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1109 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1110 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1112 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1114 if (do_ene || do_log)
1121 /* these CGLO_ options remain the same throughout the iteration */
1122 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1123 (bGStat ? CGLO_GSTAT : 0)
1126 force_flags = (GMX_FORCE_STATECHANGED |
1127 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1128 GMX_FORCE_ALLFORCES |
1130 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1131 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1132 (bDoFEP ? GMX_FORCE_DHDL : 0)
1137 if (do_per_step(step, ir->nstcalclr))
1139 force_flags |= GMX_FORCE_DO_LR;
1145 /* Now is the time to relax the shells */
1146 count = relax_shell_flexcon(fplog, cr, bVerbose, bFFscan ? step+1 : step,
1147 ir, bNS, force_flags,
1150 state, f, force_vir, mdatoms,
1151 nrnb, wcycle, graph, groups,
1152 shellfc, fr, bBornRadii, t, mu_tot,
1164 /* The coordinates (x) are shifted (to get whole molecules)
1166 * This is parallellized as well, and does communication too.
1167 * Check comments in sim_util.c
1169 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1170 state->box, state->x, &state->hist,
1171 f, force_vir, mdatoms, enerd, fcd,
1172 state->lambda, graph,
1173 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1174 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1179 mu_aver = calc_mu_aver(cr, state->x, mdatoms->chargeA,
1180 mu_tot, &top_global->mols, mdatoms, gnx, grpindex);
1183 if (bTCR && bFirstStep)
1185 tcr = init_coupling(fplog, nfile, fnm, cr, fr, mdatoms, &(top->idef));
1186 fprintf(fplog, "Done init_coupling\n");
1190 if (bVV && !bStartingFromCpt && !bRerunMD)
1191 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1193 if (ir->eI == eiVV && bInitStep)
1195 /* if using velocity verlet with full time step Ekin,
1196 * take the first half step only to compute the
1197 * virial for the first step. From there,
1198 * revert back to the initial coordinates
1199 * so that the input is actually the initial step.
1201 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1205 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1206 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1209 /* If we are using twin-range interactions where the long-range component
1210 * is only evaluated every nstcalclr>1 steps, we should do a special update
1211 * step to combine the long-range forces on these steps.
1212 * For nstcalclr=1 this is not done, since the forces would have been added
1213 * directly to the short-range forces already.
1215 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1217 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1218 f, bUpdateDoLR, fr->f_twin, fcd,
1219 ekind, M, upd, bInitStep, etrtVELOCITY1,
1220 cr, nrnb, constr, &top->idef);
1222 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1224 gmx_iterate_init(&iterate, TRUE);
1226 /* for iterations, we save these vectors, as we will be self-consistently iterating
1229 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1231 /* save the state */
1232 if (iterate.bIterationActive)
1234 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1237 bFirstIterate = TRUE;
1238 while (bFirstIterate || iterate.bIterationActive)
1240 if (iterate.bIterationActive)
1242 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1243 if (bFirstIterate && bTrotter)
1245 /* The first time through, we need a decent first estimate
1246 of veta(t+dt) to compute the constraints. Do
1247 this by computing the box volume part of the
1248 trotter integration at this time. Nothing else
1249 should be changed by this routine here. If
1250 !(first time), we start with the previous value
1253 veta_save = state->veta;
1254 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1255 vetanew = state->veta;
1256 state->veta = veta_save;
1261 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1263 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1264 state, fr->bMolPBC, graph, f,
1265 &top->idef, shake_vir,
1266 cr, nrnb, wcycle, upd, constr,
1267 TRUE, bCalcVir, vetanew);
1269 if (!bOK && !bFFscan)
1271 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1277 /* Need to unshift here if a do_force has been
1278 called in the previous step */
1279 unshift_self(graph, state->box, state->x);
1282 /* if VV, compute the pressure and constraints */
1283 /* For VV2, we strictly only need this if using pressure
1284 * control, but we really would like to have accurate pressures
1286 * Think about ways around this in the future?
1287 * For now, keep this choice in comments.
1289 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1290 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1292 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1293 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1295 bSumEkinhOld = TRUE;
1297 /* for vv, the first half of the integration actually corresponds to the previous step.
1298 So we need information from the last step in the first half of the integration */
1299 if (bGStat || do_per_step(step-1, nstglobalcomm))
1301 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1302 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1303 constr, NULL, FALSE, state->box,
1304 top_global, &pcurr, &bSumEkinhOld,
1307 | (bTemp ? CGLO_TEMPERATURE : 0)
1308 | (bPres ? CGLO_PRESSURE : 0)
1309 | (bPres ? CGLO_CONSTRAINT : 0)
1310 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1311 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1314 /* explanation of above:
1315 a) We compute Ekin at the full time step
1316 if 1) we are using the AveVel Ekin, and it's not the
1317 initial step, or 2) if we are using AveEkin, but need the full
1318 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1319 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1320 EkinAveVel because it's needed for the pressure */
1322 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1327 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1328 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1335 /* We need the kinetic energy at minus the half step for determining
1336 * the full step kinetic energy and possibly for T-coupling.*/
1337 /* This may not be quite working correctly yet . . . . */
1338 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1339 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1340 constr, NULL, FALSE, state->box,
1341 top_global, &pcurr, &bSumEkinhOld,
1342 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1347 if (iterate.bIterationActive &&
1348 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1349 state->veta, &vetanew))
1353 bFirstIterate = FALSE;
1356 if (bTrotter && !bInitStep)
1358 copy_mat(shake_vir, state->svir_prev);
1359 copy_mat(force_vir, state->fvir_prev);
1360 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1362 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1363 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1364 enerd->term[F_EKIN] = trace(ekind->ekin);
1367 /* if it's the initial step, we performed this first step just to get the constraint virial */
1368 if (bInitStep && ir->eI == eiVV)
1370 copy_rvecn(cbuf, state->v, 0, state->natoms);
1374 /* MRS -- now done iterating -- compute the conserved quantity */
1377 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1380 last_ekin = enerd->term[F_EKIN];
1382 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1384 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1386 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1389 sum_dhdl(enerd, state->lambda, ir->fepvals);
1393 /* ######## END FIRST UPDATE STEP ############## */
1394 /* ######## If doing VV, we now have v(dt) ###### */
1397 /* perform extended ensemble sampling in lambda - we don't
1398 actually move to the new state before outputting
1399 statistics, but if performing simulated tempering, we
1400 do update the velocities and the tau_t. */
1402 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, &df_history, step, mcrng, state->v, mdatoms);
1404 /* ################## START TRAJECTORY OUTPUT ################# */
1406 /* Now we have the energies and forces corresponding to the
1407 * coordinates at time t. We must output all of this before
1409 * for RerunMD t is read from input trajectory
1412 if (do_per_step(step, ir->nstxout))
1414 mdof_flags |= MDOF_X;
1416 if (do_per_step(step, ir->nstvout))
1418 mdof_flags |= MDOF_V;
1420 if (do_per_step(step, ir->nstfout))
1422 mdof_flags |= MDOF_F;
1424 if (do_per_step(step, ir->nstxtcout))
1426 mdof_flags |= MDOF_XTC;
1430 mdof_flags |= MDOF_CPT;
1434 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1437 /* Enforce writing positions and velocities at end of run */
1438 mdof_flags |= (MDOF_X | MDOF_V);
1444 fcReportProgress( ir->nsteps, step );
1447 /* sync bCPT and fc record-keeping */
1448 if (bCPT && MASTER(cr))
1450 fcRequestCheckPoint();
1454 if (mdof_flags != 0)
1456 wallcycle_start(wcycle, ewcTRAJ);
1459 if (state->flags & (1<<estLD_RNG))
1461 get_stochd_state(upd, state);
1463 if (state->flags & (1<<estMC_RNG))
1465 get_mc_state(mcrng, state);
1471 state_global->ekinstate.bUpToDate = FALSE;
1475 update_ekinstate(&state_global->ekinstate, ekind);
1476 state_global->ekinstate.bUpToDate = TRUE;
1478 update_energyhistory(&state_global->enerhist, mdebin);
1479 if (ir->efep != efepNO || ir->bSimTemp)
1481 state_global->fep_state = state->fep_state; /* MRS: seems kludgy. The code should be
1482 structured so this isn't necessary.
1483 Note this reassignment is only necessary
1484 for single threads.*/
1485 copy_df_history(&state_global->dfhist, &df_history);
1489 write_traj(fplog, cr, outf, mdof_flags, top_global,
1490 step, t, state, state_global, f, f_global, &n_xtc, &x_xtc);
1497 if (bLastStep && step_rel == ir->nsteps &&
1498 (Flags & MD_CONFOUT) && MASTER(cr) &&
1499 !bRerunMD && !bFFscan)
1501 /* x and v have been collected in write_traj,
1502 * because a checkpoint file will always be written
1505 fprintf(stderr, "\nWriting final coordinates.\n");
1508 /* Make molecules whole only for confout writing */
1509 do_pbc_mtop(fplog, ir->ePBC, state->box, top_global, state_global->x);
1511 write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm),
1512 *top_global->name, top_global,
1513 state_global->x, state_global->v,
1514 ir->ePBC, state->box);
1517 wallcycle_stop(wcycle, ewcTRAJ);
1520 /* kludge -- virial is lost with restart for NPT control. Must restart */
1521 if (bStartingFromCpt && bVV)
1523 copy_mat(state->svir_prev, shake_vir);
1524 copy_mat(state->fvir_prev, force_vir);
1526 /* ################## END TRAJECTORY OUTPUT ################ */
1528 /* Determine the wallclock run time up till now */
1529 run_time = gmx_gettime() - (double)runtime->real;
1531 /* Check whether everything is still allright */
1532 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1533 #ifdef GMX_THREAD_MPI
1538 /* this is just make gs.sig compatible with the hack
1539 of sending signals around by MPI_Reduce with together with
1541 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1543 gs.sig[eglsSTOPCOND] = 1;
1545 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1547 gs.sig[eglsSTOPCOND] = -1;
1549 /* < 0 means stop at next step, > 0 means stop at next NS step */
1553 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1554 gmx_get_signal_name(),
1555 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1559 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1560 gmx_get_signal_name(),
1561 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1563 handled_stop_condition = (int)gmx_get_stop_condition();
1565 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1566 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1567 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1569 /* Signal to terminate the run */
1570 gs.sig[eglsSTOPCOND] = 1;
1573 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1575 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1578 if (bResetCountersHalfMaxH && MASTER(cr) &&
1579 run_time > max_hours*60.0*60.0*0.495)
1581 gs.sig[eglsRESETCOUNTERS] = 1;
1584 if (ir->nstlist == -1 && !bRerunMD)
1586 /* When bGStatEveryStep=FALSE, global_stat is only called
1587 * when we check the atom displacements, not at NS steps.
1588 * This means that also the bonded interaction count check is not
1589 * performed immediately after NS. Therefore a few MD steps could
1590 * be performed with missing interactions.
1591 * But wrong energies are never written to file,
1592 * since energies are only written after global_stat
1595 if (step >= nlh.step_nscheck)
1597 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1598 nlh.scale_tot, state->x);
1602 /* This is not necessarily true,
1603 * but step_nscheck is determined quite conservatively.
1609 /* In parallel we only have to check for checkpointing in steps
1610 * where we do global communication,
1611 * otherwise the other nodes don't know.
1613 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1616 run_time >= nchkpt*cpt_period*60.0)) &&
1617 gs.set[eglsCHKPT] == 0)
1619 gs.sig[eglsCHKPT] = 1;
1622 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1627 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1629 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1631 gmx_bool bIfRandomize;
1632 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1633 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1634 if (constr && bIfRandomize)
1636 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1637 state, fr->bMolPBC, graph, f,
1638 &top->idef, tmp_vir,
1639 cr, nrnb, wcycle, upd, constr,
1640 TRUE, bCalcVir, vetanew);
1645 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1647 gmx_iterate_init(&iterate, TRUE);
1648 /* for iterations, we save these vectors, as we will be redoing the calculations */
1649 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1652 bFirstIterate = TRUE;
1653 while (bFirstIterate || iterate.bIterationActive)
1655 /* We now restore these vectors to redo the calculation with improved extended variables */
1656 if (iterate.bIterationActive)
1658 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1661 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1662 so scroll down for that logic */
1664 /* ######### START SECOND UPDATE STEP ################# */
1665 /* Box is changed in update() when we do pressure coupling,
1666 * but we should still use the old box for energy corrections and when
1667 * writing it to the energy file, so it matches the trajectory files for
1668 * the same timestep above. Make a copy in a separate array.
1670 copy_mat(state->box, lastbox);
1675 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1677 wallcycle_start(wcycle, ewcUPDATE);
1678 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1681 if (iterate.bIterationActive)
1689 /* we use a new value of scalevir to converge the iterations faster */
1690 scalevir = tracevir/trace(shake_vir);
1692 msmul(shake_vir, scalevir, shake_vir);
1693 m_add(force_vir, shake_vir, total_vir);
1694 clear_mat(shake_vir);
1696 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1697 /* We can only do Berendsen coupling after we have summed
1698 * the kinetic energy or virial. Since the happens
1699 * in global_state after update, we should only do it at
1700 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1705 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1706 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1711 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1713 /* velocity half-step update */
1714 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1715 bUpdateDoLR, fr->f_twin, fcd,
1716 ekind, M, upd, FALSE, etrtVELOCITY2,
1717 cr, nrnb, constr, &top->idef);
1720 /* Above, initialize just copies ekinh into ekin,
1721 * it doesn't copy position (for VV),
1722 * and entire integrator for MD.
1725 if (ir->eI == eiVVAK)
1727 copy_rvecn(state->x, cbuf, 0, state->natoms);
1729 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1731 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1732 bUpdateDoLR, fr->f_twin, fcd,
1733 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1734 wallcycle_stop(wcycle, ewcUPDATE);
1736 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1737 fr->bMolPBC, graph, f,
1738 &top->idef, shake_vir,
1739 cr, nrnb, wcycle, upd, constr,
1740 FALSE, bCalcVir, state->veta);
1742 if (ir->eI == eiVVAK)
1744 /* erase F_EKIN and F_TEMP here? */
1745 /* just compute the kinetic energy at the half step to perform a trotter step */
1746 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1747 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1748 constr, NULL, FALSE, lastbox,
1749 top_global, &pcurr, &bSumEkinhOld,
1750 cglo_flags | CGLO_TEMPERATURE
1752 wallcycle_start(wcycle, ewcUPDATE);
1753 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1754 /* now we know the scaling, we can compute the positions again again */
1755 copy_rvecn(cbuf, state->x, 0, state->natoms);
1757 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1759 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1760 bUpdateDoLR, fr->f_twin, fcd,
1761 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1762 wallcycle_stop(wcycle, ewcUPDATE);
1764 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1765 /* are the small terms in the shake_vir here due
1766 * to numerical errors, or are they important
1767 * physically? I'm thinking they are just errors, but not completely sure.
1768 * For now, will call without actually constraining, constr=NULL*/
1769 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1770 state, fr->bMolPBC, graph, f,
1771 &top->idef, tmp_vir,
1772 cr, nrnb, wcycle, upd, NULL,
1776 if (!bOK && !bFFscan)
1778 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1781 if (fr->bSepDVDL && fplog && do_log)
1783 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1787 /* this factor or 2 correction is necessary
1788 because half of the constraint force is removed
1789 in the vv step, so we have to double it. See
1790 the Redmine issue #1255. It is not yet clear
1791 if the factor of 2 is exact, or just a very
1792 good approximation, and this will be
1793 investigated. The next step is to see if this
1794 can be done adding a dhdl contribution from the
1795 rattle step, but this is somewhat more
1796 complicated with the current code. Will be
1797 investigated, hopefully for 4.6.3. However,
1798 this current solution is much better than
1799 having it completely wrong.
1801 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1805 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1810 /* Need to unshift here */
1811 unshift_self(graph, state->box, state->x);
1816 wallcycle_start(wcycle, ewcVSITECONSTR);
1819 shift_self(graph, state->box, state->x);
1821 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1822 top->idef.iparams, top->idef.il,
1823 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1827 unshift_self(graph, state->box, state->x);
1829 wallcycle_stop(wcycle, ewcVSITECONSTR);
1832 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1833 /* With Leap-Frog we can skip compute_globals at
1834 * non-communication steps, but we need to calculate
1835 * the kinetic energy one step before communication.
1837 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1839 if (ir->nstlist == -1 && bFirstIterate)
1841 gs.sig[eglsNABNSB] = nlh.nabnsb;
1843 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1844 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1846 bFirstIterate ? &gs : NULL,
1847 (step_rel % gs.nstms == 0) &&
1848 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1850 top_global, &pcurr, &bSumEkinhOld,
1852 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1853 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1854 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1855 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1856 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1857 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1860 if (ir->nstlist == -1 && bFirstIterate)
1862 nlh.nabnsb = gs.set[eglsNABNSB];
1863 gs.set[eglsNABNSB] = 0;
1866 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1867 /* ############# END CALC EKIN AND PRESSURE ################# */
1869 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1870 the virial that should probably be addressed eventually. state->veta has better properies,
1871 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1872 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1874 if (iterate.bIterationActive &&
1875 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1876 trace(shake_vir), &tracevir))
1880 bFirstIterate = FALSE;
1883 if (!bVV || bRerunMD)
1885 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1886 sum_dhdl(enerd, state->lambda, ir->fepvals);
1888 update_box(fplog, step, ir, mdatoms, state, f,
1889 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1891 /* ################# END UPDATE STEP 2 ################# */
1892 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1894 /* The coordinates (x) were unshifted in update */
1895 if (bFFscan && (shellfc == NULL || bConverged))
1897 if (print_forcefield(fplog, enerd->term, mdatoms->homenr,
1899 &(top_global->mols), mdatoms->massT, pres))
1903 fprintf(stderr, "\n");
1909 /* We will not sum ekinh_old,
1910 * so signal that we still have to do it.
1912 bSumEkinhOld = TRUE;
1917 /* Only do GCT when the relaxation of shells (minimization) has converged,
1918 * otherwise we might be coupling to bogus energies.
1919 * In parallel we must always do this, because the other sims might
1923 /* Since this is called with the new coordinates state->x, I assume
1924 * we want the new box state->box too. / EL 20040121
1926 do_coupling(fplog, oenv, nfile, fnm, tcr, t, step, enerd->term, fr,
1928 mdatoms, &(top->idef), mu_aver,
1929 top_global->mols.nr, cr,
1930 state->box, total_vir, pres,
1931 mu_tot, state->x, f, bConverged);
1935 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1937 /* use the directly determined last velocity, not actually the averaged half steps */
1938 if (bTrotter && ir->eI == eiVV)
1940 enerd->term[F_EKIN] = last_ekin;
1942 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1946 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1950 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1952 /* Check for excessively large energies */
1956 real etot_max = 1e200;
1958 real etot_max = 1e30;
1960 if (fabs(enerd->term[F_ETOT]) > etot_max)
1962 fprintf(stderr, "Energy too large (%g), giving up\n",
1963 enerd->term[F_ETOT]);
1966 /* ######### END PREPARING EDR OUTPUT ########### */
1968 /* Time for performance */
1969 if (((step % stepout) == 0) || bLastStep)
1971 runtime_upd_proc(runtime);
1977 gmx_bool do_dr, do_or;
1979 if (fplog && do_log && bDoExpanded)
1981 /* only needed if doing expanded ensemble */
1982 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1983 &df_history, state->fep_state, ir->nstlog, step);
1985 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1989 upd_mdebin(mdebin, bDoDHDL, TRUE,
1990 t, mdatoms->tmass, enerd, state,
1991 ir->fepvals, ir->expandedvals, lastbox,
1992 shake_vir, force_vir, total_vir, pres,
1993 ekind, mu_tot, constr);
1997 upd_mdebin_step(mdebin);
2000 do_dr = do_per_step(step, ir->nstdisreout);
2001 do_or = do_per_step(step, ir->nstorireout);
2003 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
2005 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
2007 if (ir->ePull != epullNO)
2009 pull_print_output(ir->pull, step, t);
2012 if (do_per_step(step, ir->nstlog))
2014 if (fflush(fplog) != 0)
2016 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
2022 /* Have to do this part after outputting the logfile and the edr file */
2023 state->fep_state = lamnew;
2024 for (i = 0; i < efptNR; i++)
2026 state_global->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
2029 /* Remaining runtime */
2030 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
2034 fprintf(stderr, "\n");
2036 print_time(stderr, runtime, step, ir, cr);
2039 /* Replica exchange */
2041 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
2042 do_per_step(step, repl_ex_nst))
2044 bExchanged = replica_exchange(fplog, cr, repl_ex,
2045 state_global, enerd,
2048 if (bExchanged && DOMAINDECOMP(cr))
2050 dd_partition_system(fplog, step, cr, TRUE, 1,
2051 state_global, top_global, ir,
2052 state, &f, mdatoms, top, fr,
2053 vsite, shellfc, constr,
2054 nrnb, wcycle, FALSE);
2060 bStartingFromCpt = FALSE;
2062 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
2063 /* With all integrators, except VV, we need to retain the pressure
2064 * at the current step for coupling at the next step.
2066 if ((state->flags & (1<<estPRES_PREV)) &&
2068 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
2070 /* Store the pressure in t_state for pressure coupling
2071 * at the next MD step.
2073 copy_mat(pres, state->pres_prev);
2076 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2078 if ( (membed != NULL) && (!bLastStep) )
2080 rescale_membed(step_rel, membed, state_global->x);
2087 /* read next frame from input trajectory */
2088 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
2093 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
2097 if (!bRerunMD || !rerun_fr.bStep)
2099 /* increase the MD step number */
2104 cycles = wallcycle_stop(wcycle, ewcSTEP);
2105 if (DOMAINDECOMP(cr) && wcycle)
2107 dd_cycles_add(cr->dd, cycles, ddCyclStep);
2110 if (bPMETuneRunning || bPMETuneTry)
2112 /* PME grid + cut-off optimization with GPUs or PME nodes */
2114 /* Count the total cycles over the last steps */
2115 cycles_pmes += cycles;
2117 /* We can only switch cut-off at NS steps */
2118 if (step % ir->nstlist == 0)
2120 /* PME grid + cut-off optimization with GPUs or PME nodes */
2123 if (DDMASTER(cr->dd))
2125 /* PME node load is too high, start tuning */
2126 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2128 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
2130 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2132 bPMETuneTry = FALSE;
2135 if (bPMETuneRunning)
2137 /* init_step might not be a multiple of nstlist,
2138 * but the first cycle is always skipped anyhow.
2141 pme_load_balance(pme_loadbal, cr,
2142 (bVerbose && MASTER(cr)) ? stderr : NULL,
2144 ir, state, cycles_pmes,
2145 fr->ic, fr->nbv, &fr->pmedata,
2148 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2149 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2150 fr->rlist = fr->ic->rlist;
2151 fr->rlistlong = fr->ic->rlistlong;
2152 fr->rcoulomb = fr->ic->rcoulomb;
2153 fr->rvdw = fr->ic->rvdw;
2159 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2160 gs.set[eglsRESETCOUNTERS] != 0)
2162 /* Reset all the counters related to performance over the run */
2163 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, runtime,
2164 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2165 wcycle_set_reset_counters(wcycle, -1);
2166 if (!(cr->duty & DUTY_PME))
2168 /* Tell our PME node to reset its counters */
2169 gmx_pme_send_resetcounters(cr, step);
2171 /* Correct max_hours for the elapsed time */
2172 max_hours -= run_time/(60.0*60.0);
2173 bResetCountersHalfMaxH = FALSE;
2174 gs.set[eglsRESETCOUNTERS] = 0;
2178 /* End of main MD loop */
2182 runtime_end(runtime);
2184 if (bRerunMD && MASTER(cr))
2189 if (!(cr->duty & DUTY_PME))
2191 /* Tell the PME only node to finish */
2192 gmx_pme_send_finish(cr);
2197 if (ir->nstcalcenergy > 0 && !bRerunMD)
2199 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
2200 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
2208 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2210 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)));
2211 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
2214 if (pme_loadbal != NULL)
2216 pme_loadbal_done(pme_loadbal, cr, fplog,
2217 fr->nbv != NULL && fr->nbv->bUseGPU);
2220 if (shellfc && fplog)
2222 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
2223 (nconverged*100.0)/step_rel);
2224 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
2228 if (repl_ex_nst > 0 && MASTER(cr))
2230 print_replica_exchange_statistics(fplog, repl_ex);
2233 runtime->nsteps_done = step_rel;