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;
188 gmx_groups_t *groups;
189 gmx_ekindata_t *ekind, *ekind_save;
190 gmx_shellfc_t shellfc;
191 int count, nconverged = 0;
194 gmx_bool bIonize = FALSE;
195 gmx_bool bTCR = FALSE, bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
197 gmx_bool bResetCountersHalfMaxH = FALSE;
198 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
199 gmx_bool bUpdateDoLR;
200 real mu_aver = 0, dvdl_constr;
201 int a0, a1, gnx = 0, ii;
202 atom_id *grpindex = NULL;
204 t_coupl_rec *tcr = NULL;
205 rvec *xcopy = NULL, *vcopy = NULL, *cbuf = NULL;
206 matrix boxcopy = {{0}}, lastbox;
208 real fom, oldfom, veta_save, pcurr, scalevir, tracevir;
215 real saved_conserved_quantity = 0;
220 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
221 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
222 gmx_iterate_t iterate;
223 gmx_large_int_t multisim_nsteps = -1; /* number of steps to do before first multisim
224 simulation stops. If equal to zero, don't
225 communicate any more between multisims.*/
226 /* PME load balancing data for GPU kernels */
227 pme_load_balancing_t pme_loadbal = NULL;
229 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
232 /* Temporary addition for FAHCORE checkpointing */
236 /* Check for special mdrun options */
237 bRerunMD = (Flags & MD_RERUN);
238 bIonize = (Flags & MD_IONIZE);
239 bAppend = (Flags & MD_APPENDFILES);
240 if (Flags & MD_RESETCOUNTERSHALFWAY)
244 /* Signal to reset the counters half the simulation steps. */
245 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
247 /* Signal to reset the counters halfway the simulation time. */
248 bResetCountersHalfMaxH = (max_hours > 0);
251 /* md-vv uses averaged full step velocities for T-control
252 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
253 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
255 if (bVV) /* to store the initial velocities while computing virial */
257 snew(cbuf, top_global->natoms);
259 /* all the iteratative cases - only if there are constraints */
260 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
261 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
262 false in this step. The correct value, true or false,
263 is set at each step, as it depends on the frequency of temperature
264 and pressure control.*/
265 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
269 /* Since we don't know if the frames read are related in any way,
270 * rebuild the neighborlist at every step.
273 ir->nstcalcenergy = 1;
277 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
279 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
280 bGStatEveryStep = (nstglobalcomm == 1);
282 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
285 "To reduce the energy communication with nstlist = -1\n"
286 "the neighbor list validity should not be checked at every step,\n"
287 "this means that exact integration is not guaranteed.\n"
288 "The neighbor list validity is checked after:\n"
289 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
290 "In most cases this will result in exact integration.\n"
291 "This reduces the energy communication by a factor of 2 to 3.\n"
292 "If you want less energy communication, set nstlist > 3.\n\n");
299 groups = &top_global->groups;
302 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
303 &(state_global->fep_state), lam0,
304 nrnb, top_global, &upd,
305 nfile, fnm, &outf, &mdebin,
306 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
308 clear_mat(total_vir);
310 /* Energy terms and groups */
312 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
314 if (DOMAINDECOMP(cr))
320 snew(f, top_global->natoms);
323 /* lambda Monte carlo random number generator */
326 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
328 /* copy the state into df_history */
329 copy_df_history(&df_history, &state_global->dfhist);
331 /* Kinetic energy data */
333 init_ekindata(fplog, top_global, &(ir->opts), ekind);
334 /* needed for iteration of constraints */
336 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
337 /* Copy the cos acceleration to the groups struct */
338 ekind->cosacc.cos_accel = ir->cos_accel;
340 gstat = global_stat_init(ir);
343 /* Check for polarizable models and flexible constraints */
344 shellfc = init_shell_flexcon(fplog,
345 top_global, n_flexible_constraints(constr),
346 (ir->bContinuation ||
347 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
348 NULL : state_global->x);
352 #ifdef GMX_THREAD_MPI
353 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
355 set_deform_reference_box(upd,
356 deform_init_init_step_tpx,
357 deform_init_box_tpx);
358 #ifdef GMX_THREAD_MPI
359 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
364 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
365 if ((io > 2000) && MASTER(cr))
368 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
373 if (DOMAINDECOMP(cr))
375 top = dd_init_local_top(top_global);
378 dd_init_local_state(cr->dd, state_global, state);
380 if (DDMASTER(cr->dd) && ir->nstfout)
382 snew(f_global, state_global->natoms);
389 /* Initialize the particle decomposition and split the topology */
390 top = split_system(fplog, top_global, ir, cr);
392 pd_cg_range(cr, &fr->cg0, &fr->hcg);
393 pd_at_range(cr, &a0, &a1);
397 top = gmx_mtop_generate_local_top(top_global, ir);
400 a1 = top_global->natoms;
403 forcerec_set_excl_load(fr, top, cr);
405 state = partdec_init_local_state(cr, state_global);
408 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
412 set_vsite_top(vsite, top, mdatoms, cr);
415 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
417 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
422 make_local_shells(cr, mdatoms, shellfc);
425 setup_bonded_threading(fr, &top->idef);
427 if (ir->pull && PAR(cr))
429 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
433 if (DOMAINDECOMP(cr))
435 /* Distribute the charge groups over the nodes from the master node */
436 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
437 state_global, top_global, ir,
438 state, &f, mdatoms, top, fr,
439 vsite, shellfc, constr,
440 nrnb, wcycle, FALSE);
444 update_mdatoms(mdatoms, state->lambda[efptMASS]);
446 if (opt2bSet("-cpi", nfile, fnm))
448 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
452 bStateFromCP = FALSE;
459 /* Update mdebin with energy history if appending to output files */
460 if (Flags & MD_APPENDFILES)
462 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
466 /* We might have read an energy history from checkpoint,
467 * free the allocated memory and reset the counts.
469 done_energyhistory(&state_global->enerhist);
470 init_energyhistory(&state_global->enerhist);
473 /* Set the initial energy history in state by updating once */
474 update_energyhistory(&state_global->enerhist, mdebin);
477 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
479 /* Set the random state if we read a checkpoint file */
480 set_stochd_state(upd, state);
483 if (state->flags & (1<<estMC_RNG))
485 set_mc_state(mcrng, state);
488 /* Initialize constraints */
491 if (!DOMAINDECOMP(cr))
493 set_constraints(constr, top, ir, mdatoms, cr);
497 /* Check whether we have to GCT stuff */
498 bTCR = ftp2bSet(efGCT, nfile, fnm);
503 fprintf(stderr, "Will do General Coupling Theory!\n");
505 gnx = top_global->mols.nr;
507 for (i = 0; (i < gnx); i++)
515 /* We need to be sure replica exchange can only occur
516 * when the energies are current */
517 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
518 "repl_ex_nst", &repl_ex_nst);
519 /* This check needs to happen before inter-simulation
520 * signals are initialized, too */
522 if (repl_ex_nst > 0 && MASTER(cr))
524 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
525 repl_ex_nst, repl_ex_nex, repl_ex_seed);
528 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
529 * With perturbed charges with soft-core we should not change the cut-off.
531 if ((Flags & MD_TUNEPME) &&
532 EEL_PME(fr->eeltype) &&
533 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
534 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
537 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
539 if (cr->duty & DUTY_PME)
541 /* Start tuning right away, as we can't measure the load */
542 bPMETuneRunning = TRUE;
546 /* Separate PME nodes, we can measure the PP/PME load balance */
551 if (!ir->bContinuation && !bRerunMD)
553 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
555 /* Set the velocities of frozen particles to zero */
556 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
558 for (m = 0; m < DIM; m++)
560 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
570 /* Constrain the initial coordinates and velocities */
571 do_constrain_first(fplog, constr, ir, mdatoms, state,
576 /* Construct the virtual sites for the initial configuration */
577 construct_vsites(vsite, state->x, ir->delta_t, NULL,
578 top->idef.iparams, top->idef.il,
579 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
585 /* set free energy calculation frequency as the minimum
586 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
587 nstfep = ir->fepvals->nstdhdl;
590 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
594 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
597 /* I'm assuming we need global communication the first time! MRS */
598 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
599 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
600 | (bVV ? CGLO_PRESSURE : 0)
601 | (bVV ? CGLO_CONSTRAINT : 0)
602 | (bRerunMD ? CGLO_RERUNMD : 0)
603 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
605 bSumEkinhOld = FALSE;
606 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
607 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
608 constr, NULL, FALSE, state->box,
609 top_global, &pcurr, &bSumEkinhOld, cglo_flags);
610 if (ir->eI == eiVVAK)
612 /* a second call to get the half step temperature initialized as well */
613 /* we do the same call as above, but turn the pressure off -- internally to
614 compute_globals, this is recognized as a velocity verlet half-step
615 kinetic energy calculation. This minimized excess variables, but
616 perhaps loses some logic?*/
618 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
619 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
620 constr, NULL, FALSE, state->box,
621 top_global, &pcurr, &bSumEkinhOld,
622 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
625 /* Calculate the initial half step temperature, and save the ekinh_old */
626 if (!(Flags & MD_STARTFROMCPT))
628 for (i = 0; (i < ir->opts.ngtc); i++)
630 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
635 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
636 and there is no previous step */
639 /* if using an iterative algorithm, we need to create a working directory for the state. */
642 bufstate = init_bufstate(state);
645 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
646 temperature control */
647 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
651 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
654 "RMS relative constraint deviation after constraining: %.2e\n",
655 constr_rmsd(constr, FALSE));
657 if (EI_STATE_VELOCITY(ir->eI))
659 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
663 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
664 " input trajectory '%s'\n\n",
665 *(top_global->name), opt2fn("-rerun", nfile, fnm));
668 fprintf(stderr, "Calculated time to finish depends on nsteps from "
669 "run input file,\nwhich may not correspond to the time "
670 "needed to process input trajectory.\n\n");
676 fprintf(stderr, "starting mdrun '%s'\n",
677 *(top_global->name));
680 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
684 sprintf(tbuf, "%s", "infinite");
686 if (ir->init_step > 0)
688 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
689 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
690 gmx_step_str(ir->init_step, sbuf2),
691 ir->init_step*ir->delta_t);
695 fprintf(stderr, "%s steps, %s ps.\n",
696 gmx_step_str(ir->nsteps, sbuf), tbuf);
699 fprintf(fplog, "\n");
702 /* Set and write start time */
703 runtime_start(runtime);
704 print_date_and_time(fplog, cr->nodeid, "Started mdrun", runtime);
705 wallcycle_start(wcycle, ewcRUN);
708 fprintf(fplog, "\n");
711 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
713 chkpt_ret = fcCheckPointParallel( cr->nodeid,
717 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
722 /***********************************************************
726 ************************************************************/
728 /* if rerunMD then read coordinates and velocities from input trajectory */
731 if (getenv("GMX_FORCE_UPDATE"))
739 bNotLastFrame = read_first_frame(oenv, &status,
740 opt2fn("-rerun", nfile, fnm),
741 &rerun_fr, TRX_NEED_X | TRX_READ_V);
742 if (rerun_fr.natoms != top_global->natoms)
745 "Number of atoms in trajectory (%d) does not match the "
746 "run input file (%d)\n",
747 rerun_fr.natoms, top_global->natoms);
749 if (ir->ePBC != epbcNONE)
753 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);
755 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
757 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
764 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
767 if (ir->ePBC != epbcNONE)
769 /* Set the shift vectors.
770 * Necessary here when have a static box different from the tpr box.
772 calc_shifts(rerun_fr.box, fr->shift_vec);
776 /* loop over MD steps or if rerunMD to end of input trajectory */
778 /* Skip the first Nose-Hoover integration when we get the state from tpx */
779 bStateFromTPX = !bStateFromCP;
780 bInitStep = bFirstStep && (bStateFromTPX || bVV);
781 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
783 bSumEkinhOld = FALSE;
786 init_global_signals(&gs, cr, ir, repl_ex_nst);
788 step = ir->init_step;
791 if (ir->nstlist == -1)
793 init_nlistheuristics(&nlh, bGStatEveryStep, step);
796 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
798 /* check how many steps are left in other sims */
799 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
803 /* and stop now if we should */
804 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
805 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
806 while (!bLastStep || (bRerunMD && bNotLastFrame))
809 wallcycle_start(wcycle, ewcSTEP);
815 step = rerun_fr.step;
816 step_rel = step - ir->init_step;
829 bLastStep = (step_rel == ir->nsteps);
830 t = t0 + step*ir->delta_t;
833 if (ir->efep != efepNO || ir->bSimTemp)
835 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
836 requiring different logic. */
838 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
839 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
840 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
841 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
846 update_annealing_target_temp(&(ir->opts), t);
851 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
853 for (i = 0; i < state_global->natoms; i++)
855 copy_rvec(rerun_fr.x[i], state_global->x[i]);
859 for (i = 0; i < state_global->natoms; i++)
861 copy_rvec(rerun_fr.v[i], state_global->v[i]);
866 for (i = 0; i < state_global->natoms; i++)
868 clear_rvec(state_global->v[i]);
872 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
873 " Ekin, temperature and pressure are incorrect,\n"
874 " the virial will be incorrect when constraints are present.\n"
876 bRerunWarnNoV = FALSE;
880 copy_mat(rerun_fr.box, state_global->box);
881 copy_mat(state_global->box, state->box);
883 if (vsite && (Flags & MD_RERUN_VSITE))
885 if (DOMAINDECOMP(cr))
887 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
891 /* Following is necessary because the graph may get out of sync
892 * with the coordinates if we only have every N'th coordinate set
894 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
895 shift_self(graph, state->box, state->x);
897 construct_vsites(vsite, state->x, ir->delta_t, state->v,
898 top->idef.iparams, top->idef.il,
899 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
902 unshift_self(graph, state->box, state->x);
907 /* Stop Center of Mass motion */
908 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
912 /* for rerun MD always do Neighbour Searching */
913 bNS = (bFirstStep || ir->nstlist != 0);
918 /* Determine whether or not to do Neighbour Searching and LR */
919 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
921 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
922 (ir->nstlist == -1 && nlh.nabnsb > 0));
924 if (bNS && ir->nstlist == -1)
926 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
930 /* check whether we should stop because another simulation has
934 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
935 (multisim_nsteps != ir->nsteps) )
942 "Stopping simulation %d because another one has finished\n",
946 gs.sig[eglsCHKPT] = 1;
951 /* < 0 means stop at next step, > 0 means stop at next NS step */
952 if ( (gs.set[eglsSTOPCOND] < 0) ||
953 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
958 /* Determine whether or not to update the Born radii if doing GB */
959 bBornRadii = bFirstStep;
960 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
965 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
966 do_verbose = bVerbose &&
967 (step % stepout == 0 || bFirstStep || bLastStep);
969 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
977 bMasterState = FALSE;
978 /* Correct the new box if it is too skewed */
979 if (DYNAMIC_BOX(*ir))
981 if (correct_box(fplog, step, state->box, graph))
986 if (DOMAINDECOMP(cr) && bMasterState)
988 dd_collect_state(cr->dd, state, state_global);
992 if (DOMAINDECOMP(cr))
994 /* Repartition the domain decomposition */
995 wallcycle_start(wcycle, ewcDOMDEC);
996 dd_partition_system(fplog, step, cr,
997 bMasterState, nstglobalcomm,
998 state_global, top_global, ir,
999 state, &f, mdatoms, top, fr,
1000 vsite, shellfc, constr,
1002 do_verbose && !bPMETuneRunning);
1003 wallcycle_stop(wcycle, ewcDOMDEC);
1004 /* If using an iterative integrator, reallocate space to match the decomposition */
1008 if (MASTER(cr) && do_log)
1010 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
1013 if (ir->efep != efepNO)
1015 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1018 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1021 /* We need the kinetic energy at minus the half step for determining
1022 * the full step kinetic energy and possibly for T-coupling.*/
1023 /* This may not be quite working correctly yet . . . . */
1024 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1025 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1026 constr, NULL, FALSE, state->box,
1027 top_global, &pcurr, &bSumEkinhOld,
1028 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1030 clear_mat(force_vir);
1032 /* Ionize the atoms if necessary */
1035 ionize(fplog, oenv, mdatoms, top_global, t, ir, state->x, state->v,
1036 mdatoms->start, mdatoms->start+mdatoms->homenr, state->box, cr);
1039 /* We write a checkpoint at this MD step when:
1040 * either at an NS step when we signalled through gs,
1041 * or at the last step (but not when we do not want confout),
1042 * but never at the first step or with rerun.
1044 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1045 (bLastStep && (Flags & MD_CONFOUT))) &&
1046 step > ir->init_step && !bRerunMD);
1049 gs.set[eglsCHKPT] = 0;
1052 /* Determine the energy and pressure:
1053 * at nstcalcenergy steps and at energy output steps (set below).
1055 if (EI_VV(ir->eI) && (!bInitStep))
1057 /* for vv, the first half of the integration actually corresponds
1058 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1059 but the virial needs to be calculated on both the current step and the 'next' step. Future
1060 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1062 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1063 bCalcVir = bCalcEner ||
1064 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1068 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1069 bCalcVir = bCalcEner ||
1070 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1073 /* Do we need global communication ? */
1074 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1075 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1076 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1078 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1080 if (do_ene || do_log)
1087 /* these CGLO_ options remain the same throughout the iteration */
1088 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1089 (bGStat ? CGLO_GSTAT : 0)
1092 force_flags = (GMX_FORCE_STATECHANGED |
1093 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1094 GMX_FORCE_ALLFORCES |
1096 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1097 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1098 (bDoFEP ? GMX_FORCE_DHDL : 0)
1103 if (do_per_step(step, ir->nstcalclr))
1105 force_flags |= GMX_FORCE_DO_LR;
1111 /* Now is the time to relax the shells */
1112 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1113 ir, bNS, force_flags,
1116 state, f, force_vir, mdatoms,
1117 nrnb, wcycle, graph, groups,
1118 shellfc, fr, bBornRadii, t, mu_tot,
1130 /* The coordinates (x) are shifted (to get whole molecules)
1132 * This is parallellized as well, and does communication too.
1133 * Check comments in sim_util.c
1135 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1136 state->box, state->x, &state->hist,
1137 f, force_vir, mdatoms, enerd, fcd,
1138 state->lambda, graph,
1139 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1140 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1145 mu_aver = calc_mu_aver(cr, state->x, mdatoms->chargeA,
1146 mu_tot, &top_global->mols, mdatoms, gnx, grpindex);
1149 if (bTCR && bFirstStep)
1151 tcr = init_coupling(fplog, nfile, fnm, cr, fr, mdatoms, &(top->idef));
1152 fprintf(fplog, "Done init_coupling\n");
1156 if (bVV && !bStartingFromCpt && !bRerunMD)
1157 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1159 if (ir->eI == eiVV && bInitStep)
1161 /* if using velocity verlet with full time step Ekin,
1162 * take the first half step only to compute the
1163 * virial for the first step. From there,
1164 * revert back to the initial coordinates
1165 * so that the input is actually the initial step.
1167 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1171 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1172 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1175 /* If we are using twin-range interactions where the long-range component
1176 * is only evaluated every nstcalclr>1 steps, we should do a special update
1177 * step to combine the long-range forces on these steps.
1178 * For nstcalclr=1 this is not done, since the forces would have been added
1179 * directly to the short-range forces already.
1181 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1183 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1184 f, bUpdateDoLR, fr->f_twin, fcd,
1185 ekind, M, upd, bInitStep, etrtVELOCITY1,
1186 cr, nrnb, constr, &top->idef);
1188 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1190 gmx_iterate_init(&iterate, TRUE);
1192 /* for iterations, we save these vectors, as we will be self-consistently iterating
1195 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1197 /* save the state */
1198 if (iterate.bIterationActive)
1200 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1203 bFirstIterate = TRUE;
1204 while (bFirstIterate || iterate.bIterationActive)
1206 if (iterate.bIterationActive)
1208 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1209 if (bFirstIterate && bTrotter)
1211 /* The first time through, we need a decent first estimate
1212 of veta(t+dt) to compute the constraints. Do
1213 this by computing the box volume part of the
1214 trotter integration at this time. Nothing else
1215 should be changed by this routine here. If
1216 !(first time), we start with the previous value
1219 veta_save = state->veta;
1220 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1221 vetanew = state->veta;
1222 state->veta = veta_save;
1227 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1229 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1230 state, fr->bMolPBC, graph, f,
1231 &top->idef, shake_vir,
1232 cr, nrnb, wcycle, upd, constr,
1233 TRUE, bCalcVir, vetanew);
1237 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1243 /* Need to unshift here if a do_force has been
1244 called in the previous step */
1245 unshift_self(graph, state->box, state->x);
1248 /* if VV, compute the pressure and constraints */
1249 /* For VV2, we strictly only need this if using pressure
1250 * control, but we really would like to have accurate pressures
1252 * Think about ways around this in the future?
1253 * For now, keep this choice in comments.
1255 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1256 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1258 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1259 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1261 bSumEkinhOld = TRUE;
1263 /* for vv, the first half of the integration actually corresponds to the previous step.
1264 So we need information from the last step in the first half of the integration */
1265 if (bGStat || do_per_step(step-1, nstglobalcomm))
1267 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1268 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1269 constr, NULL, FALSE, state->box,
1270 top_global, &pcurr, &bSumEkinhOld,
1273 | (bTemp ? CGLO_TEMPERATURE : 0)
1274 | (bPres ? CGLO_PRESSURE : 0)
1275 | (bPres ? CGLO_CONSTRAINT : 0)
1276 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1277 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1280 /* explanation of above:
1281 a) We compute Ekin at the full time step
1282 if 1) we are using the AveVel Ekin, and it's not the
1283 initial step, or 2) if we are using AveEkin, but need the full
1284 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1285 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1286 EkinAveVel because it's needed for the pressure */
1288 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1293 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1294 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1301 /* We need the kinetic energy at minus the half step for determining
1302 * the full step kinetic energy and possibly for T-coupling.*/
1303 /* This may not be quite working correctly yet . . . . */
1304 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1305 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1306 constr, NULL, FALSE, state->box,
1307 top_global, &pcurr, &bSumEkinhOld,
1308 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1313 if (iterate.bIterationActive &&
1314 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1315 state->veta, &vetanew))
1319 bFirstIterate = FALSE;
1322 if (bTrotter && !bInitStep)
1324 copy_mat(shake_vir, state->svir_prev);
1325 copy_mat(force_vir, state->fvir_prev);
1326 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1328 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1329 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1330 enerd->term[F_EKIN] = trace(ekind->ekin);
1333 /* if it's the initial step, we performed this first step just to get the constraint virial */
1334 if (bInitStep && ir->eI == eiVV)
1336 copy_rvecn(cbuf, state->v, 0, state->natoms);
1340 /* MRS -- now done iterating -- compute the conserved quantity */
1343 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1346 last_ekin = enerd->term[F_EKIN];
1348 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1350 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1352 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1355 sum_dhdl(enerd, state->lambda, ir->fepvals);
1359 /* ######## END FIRST UPDATE STEP ############## */
1360 /* ######## If doing VV, we now have v(dt) ###### */
1363 /* perform extended ensemble sampling in lambda - we don't
1364 actually move to the new state before outputting
1365 statistics, but if performing simulated tempering, we
1366 do update the velocities and the tau_t. */
1368 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, &df_history, step, mcrng, state->v, mdatoms);
1370 /* ################## START TRAJECTORY OUTPUT ################# */
1372 /* Now we have the energies and forces corresponding to the
1373 * coordinates at time t. We must output all of this before
1375 * for RerunMD t is read from input trajectory
1378 if (do_per_step(step, ir->nstxout))
1380 mdof_flags |= MDOF_X;
1382 if (do_per_step(step, ir->nstvout))
1384 mdof_flags |= MDOF_V;
1386 if (do_per_step(step, ir->nstfout))
1388 mdof_flags |= MDOF_F;
1390 if (do_per_step(step, ir->nstxtcout))
1392 mdof_flags |= MDOF_XTC;
1396 mdof_flags |= MDOF_CPT;
1400 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1403 /* Enforce writing positions and velocities at end of run */
1404 mdof_flags |= (MDOF_X | MDOF_V);
1410 fcReportProgress( ir->nsteps, step );
1413 /* sync bCPT and fc record-keeping */
1414 if (bCPT && MASTER(cr))
1416 fcRequestCheckPoint();
1420 if (mdof_flags != 0)
1422 wallcycle_start(wcycle, ewcTRAJ);
1425 if (state->flags & (1<<estLD_RNG))
1427 get_stochd_state(upd, state);
1429 if (state->flags & (1<<estMC_RNG))
1431 get_mc_state(mcrng, state);
1437 state_global->ekinstate.bUpToDate = FALSE;
1441 update_ekinstate(&state_global->ekinstate, ekind);
1442 state_global->ekinstate.bUpToDate = TRUE;
1444 update_energyhistory(&state_global->enerhist, mdebin);
1445 if (ir->efep != efepNO || ir->bSimTemp)
1447 state_global->fep_state = state->fep_state; /* MRS: seems kludgy. The code should be
1448 structured so this isn't necessary.
1449 Note this reassignment is only necessary
1450 for single threads.*/
1451 copy_df_history(&state_global->dfhist, &df_history);
1455 write_traj(fplog, cr, outf, mdof_flags, top_global,
1456 step, t, state, state_global, f, f_global, &n_xtc, &x_xtc);
1463 if (bLastStep && step_rel == ir->nsteps &&
1464 (Flags & MD_CONFOUT) && MASTER(cr) &&
1467 /* x and v have been collected in write_traj,
1468 * because a checkpoint file will always be written
1471 fprintf(stderr, "\nWriting final coordinates.\n");
1474 /* Make molecules whole only for confout writing */
1475 do_pbc_mtop(fplog, ir->ePBC, state->box, top_global, state_global->x);
1477 write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm),
1478 *top_global->name, top_global,
1479 state_global->x, state_global->v,
1480 ir->ePBC, state->box);
1483 wallcycle_stop(wcycle, ewcTRAJ);
1486 /* kludge -- virial is lost with restart for NPT control. Must restart */
1487 if (bStartingFromCpt && bVV)
1489 copy_mat(state->svir_prev, shake_vir);
1490 copy_mat(state->fvir_prev, force_vir);
1492 /* ################## END TRAJECTORY OUTPUT ################ */
1494 /* Determine the wallclock run time up till now */
1495 run_time = gmx_gettime() - (double)runtime->real;
1497 /* Check whether everything is still allright */
1498 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1499 #ifdef GMX_THREAD_MPI
1504 /* this is just make gs.sig compatible with the hack
1505 of sending signals around by MPI_Reduce with together with
1507 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1509 gs.sig[eglsSTOPCOND] = 1;
1511 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1513 gs.sig[eglsSTOPCOND] = -1;
1515 /* < 0 means stop at next step, > 0 means stop at next NS step */
1519 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1520 gmx_get_signal_name(),
1521 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1525 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1526 gmx_get_signal_name(),
1527 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1529 handled_stop_condition = (int)gmx_get_stop_condition();
1531 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1532 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1533 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1535 /* Signal to terminate the run */
1536 gs.sig[eglsSTOPCOND] = 1;
1539 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1541 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1544 if (bResetCountersHalfMaxH && MASTER(cr) &&
1545 run_time > max_hours*60.0*60.0*0.495)
1547 gs.sig[eglsRESETCOUNTERS] = 1;
1550 if (ir->nstlist == -1 && !bRerunMD)
1552 /* When bGStatEveryStep=FALSE, global_stat is only called
1553 * when we check the atom displacements, not at NS steps.
1554 * This means that also the bonded interaction count check is not
1555 * performed immediately after NS. Therefore a few MD steps could
1556 * be performed with missing interactions.
1557 * But wrong energies are never written to file,
1558 * since energies are only written after global_stat
1561 if (step >= nlh.step_nscheck)
1563 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1564 nlh.scale_tot, state->x);
1568 /* This is not necessarily true,
1569 * but step_nscheck is determined quite conservatively.
1575 /* In parallel we only have to check for checkpointing in steps
1576 * where we do global communication,
1577 * otherwise the other nodes don't know.
1579 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1582 run_time >= nchkpt*cpt_period*60.0)) &&
1583 gs.set[eglsCHKPT] == 0)
1585 gs.sig[eglsCHKPT] = 1;
1588 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1593 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1595 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1597 gmx_bool bIfRandomize;
1598 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1599 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1600 if (constr && bIfRandomize)
1602 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1603 state, fr->bMolPBC, graph, f,
1604 &top->idef, tmp_vir,
1605 cr, nrnb, wcycle, upd, constr,
1606 TRUE, bCalcVir, vetanew);
1611 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1613 gmx_iterate_init(&iterate, TRUE);
1614 /* for iterations, we save these vectors, as we will be redoing the calculations */
1615 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1618 bFirstIterate = TRUE;
1619 while (bFirstIterate || iterate.bIterationActive)
1621 /* We now restore these vectors to redo the calculation with improved extended variables */
1622 if (iterate.bIterationActive)
1624 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1627 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1628 so scroll down for that logic */
1630 /* ######### START SECOND UPDATE STEP ################# */
1631 /* Box is changed in update() when we do pressure coupling,
1632 * but we should still use the old box for energy corrections and when
1633 * writing it to the energy file, so it matches the trajectory files for
1634 * the same timestep above. Make a copy in a separate array.
1636 copy_mat(state->box, lastbox);
1641 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1643 wallcycle_start(wcycle, ewcUPDATE);
1644 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1647 if (iterate.bIterationActive)
1655 /* we use a new value of scalevir to converge the iterations faster */
1656 scalevir = tracevir/trace(shake_vir);
1658 msmul(shake_vir, scalevir, shake_vir);
1659 m_add(force_vir, shake_vir, total_vir);
1660 clear_mat(shake_vir);
1662 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1663 /* We can only do Berendsen coupling after we have summed
1664 * the kinetic energy or virial. Since the happens
1665 * in global_state after update, we should only do it at
1666 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1671 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1672 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1677 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1679 /* velocity half-step update */
1680 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1681 bUpdateDoLR, fr->f_twin, fcd,
1682 ekind, M, upd, FALSE, etrtVELOCITY2,
1683 cr, nrnb, constr, &top->idef);
1686 /* Above, initialize just copies ekinh into ekin,
1687 * it doesn't copy position (for VV),
1688 * and entire integrator for MD.
1691 if (ir->eI == eiVVAK)
1693 copy_rvecn(state->x, cbuf, 0, state->natoms);
1695 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1697 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1698 bUpdateDoLR, fr->f_twin, fcd,
1699 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1700 wallcycle_stop(wcycle, ewcUPDATE);
1702 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1703 fr->bMolPBC, graph, f,
1704 &top->idef, shake_vir,
1705 cr, nrnb, wcycle, upd, constr,
1706 FALSE, bCalcVir, state->veta);
1708 if (ir->eI == eiVVAK)
1710 /* erase F_EKIN and F_TEMP here? */
1711 /* just compute the kinetic energy at the half step to perform a trotter step */
1712 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1713 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1714 constr, NULL, FALSE, lastbox,
1715 top_global, &pcurr, &bSumEkinhOld,
1716 cglo_flags | CGLO_TEMPERATURE
1718 wallcycle_start(wcycle, ewcUPDATE);
1719 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1720 /* now we know the scaling, we can compute the positions again again */
1721 copy_rvecn(cbuf, state->x, 0, state->natoms);
1723 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1725 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1726 bUpdateDoLR, fr->f_twin, fcd,
1727 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1728 wallcycle_stop(wcycle, ewcUPDATE);
1730 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1731 /* are the small terms in the shake_vir here due
1732 * to numerical errors, or are they important
1733 * physically? I'm thinking they are just errors, but not completely sure.
1734 * For now, will call without actually constraining, constr=NULL*/
1735 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1736 state, fr->bMolPBC, graph, f,
1737 &top->idef, tmp_vir,
1738 cr, nrnb, wcycle, upd, NULL,
1744 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1747 if (fr->bSepDVDL && fplog && do_log)
1749 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1753 /* this factor or 2 correction is necessary
1754 because half of the constraint force is removed
1755 in the vv step, so we have to double it. See
1756 the Redmine issue #1255. It is not yet clear
1757 if the factor of 2 is exact, or just a very
1758 good approximation, and this will be
1759 investigated. The next step is to see if this
1760 can be done adding a dhdl contribution from the
1761 rattle step, but this is somewhat more
1762 complicated with the current code. Will be
1763 investigated, hopefully for 4.6.3. However,
1764 this current solution is much better than
1765 having it completely wrong.
1767 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1771 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1776 /* Need to unshift here */
1777 unshift_self(graph, state->box, state->x);
1782 wallcycle_start(wcycle, ewcVSITECONSTR);
1785 shift_self(graph, state->box, state->x);
1787 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1788 top->idef.iparams, top->idef.il,
1789 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1793 unshift_self(graph, state->box, state->x);
1795 wallcycle_stop(wcycle, ewcVSITECONSTR);
1798 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1799 /* With Leap-Frog we can skip compute_globals at
1800 * non-communication steps, but we need to calculate
1801 * the kinetic energy one step before communication.
1803 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1805 if (ir->nstlist == -1 && bFirstIterate)
1807 gs.sig[eglsNABNSB] = nlh.nabnsb;
1809 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1810 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1812 bFirstIterate ? &gs : NULL,
1813 (step_rel % gs.nstms == 0) &&
1814 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1816 top_global, &pcurr, &bSumEkinhOld,
1818 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1819 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1820 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1821 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1822 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1823 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1826 if (ir->nstlist == -1 && bFirstIterate)
1828 nlh.nabnsb = gs.set[eglsNABNSB];
1829 gs.set[eglsNABNSB] = 0;
1832 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1833 /* ############# END CALC EKIN AND PRESSURE ################# */
1835 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1836 the virial that should probably be addressed eventually. state->veta has better properies,
1837 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1838 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1840 if (iterate.bIterationActive &&
1841 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1842 trace(shake_vir), &tracevir))
1846 bFirstIterate = FALSE;
1849 if (!bVV || bRerunMD)
1851 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1852 sum_dhdl(enerd, state->lambda, ir->fepvals);
1854 update_box(fplog, step, ir, mdatoms, state, f,
1855 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1857 /* ################# END UPDATE STEP 2 ################# */
1858 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1860 /* The coordinates (x) were unshifted in update */
1863 /* We will not sum ekinh_old,
1864 * so signal that we still have to do it.
1866 bSumEkinhOld = TRUE;
1871 /* Only do GCT when the relaxation of shells (minimization) has converged,
1872 * otherwise we might be coupling to bogus energies.
1873 * In parallel we must always do this, because the other sims might
1877 /* Since this is called with the new coordinates state->x, I assume
1878 * we want the new box state->box too. / EL 20040121
1880 do_coupling(fplog, oenv, nfile, fnm, tcr, t, step, enerd->term, fr,
1882 mdatoms, &(top->idef), mu_aver,
1883 top_global->mols.nr, cr,
1884 state->box, total_vir, pres,
1885 mu_tot, state->x, f, bConverged);
1889 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1891 /* use the directly determined last velocity, not actually the averaged half steps */
1892 if (bTrotter && ir->eI == eiVV)
1894 enerd->term[F_EKIN] = last_ekin;
1896 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1900 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1904 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1906 /* Check for excessively large energies */
1910 real etot_max = 1e200;
1912 real etot_max = 1e30;
1914 if (fabs(enerd->term[F_ETOT]) > etot_max)
1916 fprintf(stderr, "Energy too large (%g), giving up\n",
1917 enerd->term[F_ETOT]);
1920 /* ######### END PREPARING EDR OUTPUT ########### */
1922 /* Time for performance */
1923 if (((step % stepout) == 0) || bLastStep)
1925 runtime_upd_proc(runtime);
1931 gmx_bool do_dr, do_or;
1933 if (fplog && do_log && bDoExpanded)
1935 /* only needed if doing expanded ensemble */
1936 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1937 &df_history, state->fep_state, ir->nstlog, step);
1939 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1943 upd_mdebin(mdebin, bDoDHDL, TRUE,
1944 t, mdatoms->tmass, enerd, state,
1945 ir->fepvals, ir->expandedvals, lastbox,
1946 shake_vir, force_vir, total_vir, pres,
1947 ekind, mu_tot, constr);
1951 upd_mdebin_step(mdebin);
1954 do_dr = do_per_step(step, ir->nstdisreout);
1955 do_or = do_per_step(step, ir->nstorireout);
1957 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
1959 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1961 if (ir->ePull != epullNO)
1963 pull_print_output(ir->pull, step, t);
1966 if (do_per_step(step, ir->nstlog))
1968 if (fflush(fplog) != 0)
1970 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1976 /* Have to do this part after outputting the logfile and the edr file */
1977 state->fep_state = lamnew;
1978 for (i = 0; i < efptNR; i++)
1980 state_global->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
1983 /* Remaining runtime */
1984 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1988 fprintf(stderr, "\n");
1990 print_time(stderr, runtime, step, ir, cr);
1993 /* Replica exchange */
1995 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1996 do_per_step(step, repl_ex_nst))
1998 bExchanged = replica_exchange(fplog, cr, repl_ex,
1999 state_global, enerd,
2002 if (bExchanged && DOMAINDECOMP(cr))
2004 dd_partition_system(fplog, step, cr, TRUE, 1,
2005 state_global, top_global, ir,
2006 state, &f, mdatoms, top, fr,
2007 vsite, shellfc, constr,
2008 nrnb, wcycle, FALSE);
2014 bStartingFromCpt = FALSE;
2016 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
2017 /* With all integrators, except VV, we need to retain the pressure
2018 * at the current step for coupling at the next step.
2020 if ((state->flags & (1<<estPRES_PREV)) &&
2022 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
2024 /* Store the pressure in t_state for pressure coupling
2025 * at the next MD step.
2027 copy_mat(pres, state->pres_prev);
2030 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2032 if ( (membed != NULL) && (!bLastStep) )
2034 rescale_membed(step_rel, membed, state_global->x);
2041 /* read next frame from input trajectory */
2042 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
2047 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
2051 if (!bRerunMD || !rerun_fr.bStep)
2053 /* increase the MD step number */
2058 cycles = wallcycle_stop(wcycle, ewcSTEP);
2059 if (DOMAINDECOMP(cr) && wcycle)
2061 dd_cycles_add(cr->dd, cycles, ddCyclStep);
2064 if (bPMETuneRunning || bPMETuneTry)
2066 /* PME grid + cut-off optimization with GPUs or PME nodes */
2068 /* Count the total cycles over the last steps */
2069 cycles_pmes += cycles;
2071 /* We can only switch cut-off at NS steps */
2072 if (step % ir->nstlist == 0)
2074 /* PME grid + cut-off optimization with GPUs or PME nodes */
2077 if (DDMASTER(cr->dd))
2079 /* PME node load is too high, start tuning */
2080 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2082 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
2084 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2086 bPMETuneTry = FALSE;
2089 if (bPMETuneRunning)
2091 /* init_step might not be a multiple of nstlist,
2092 * but the first cycle is always skipped anyhow.
2095 pme_load_balance(pme_loadbal, cr,
2096 (bVerbose && MASTER(cr)) ? stderr : NULL,
2098 ir, state, cycles_pmes,
2099 fr->ic, fr->nbv, &fr->pmedata,
2102 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2103 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2104 fr->rlist = fr->ic->rlist;
2105 fr->rlistlong = fr->ic->rlistlong;
2106 fr->rcoulomb = fr->ic->rcoulomb;
2107 fr->rvdw = fr->ic->rvdw;
2113 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2114 gs.set[eglsRESETCOUNTERS] != 0)
2116 /* Reset all the counters related to performance over the run */
2117 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, runtime,
2118 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2119 wcycle_set_reset_counters(wcycle, -1);
2120 if (!(cr->duty & DUTY_PME))
2122 /* Tell our PME node to reset its counters */
2123 gmx_pme_send_resetcounters(cr, step);
2125 /* Correct max_hours for the elapsed time */
2126 max_hours -= run_time/(60.0*60.0);
2127 bResetCountersHalfMaxH = FALSE;
2128 gs.set[eglsRESETCOUNTERS] = 0;
2132 /* End of main MD loop */
2136 runtime_end(runtime);
2138 if (bRerunMD && MASTER(cr))
2143 if (!(cr->duty & DUTY_PME))
2145 /* Tell the PME only node to finish */
2146 gmx_pme_send_finish(cr);
2151 if (ir->nstcalcenergy > 0 && !bRerunMD)
2153 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
2154 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
2162 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2164 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)));
2165 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
2168 if (pme_loadbal != NULL)
2170 pme_loadbal_done(pme_loadbal, cr, fplog,
2171 fr->nbv != NULL && fr->nbv->bUseGPU);
2174 if (shellfc && fplog)
2176 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
2177 (nconverged*100.0)/step_rel);
2178 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
2182 if (repl_ex_nst > 0 && MASTER(cr))
2184 print_replica_exchange_statistics(fplog, repl_ex);
2187 runtime->nsteps_done = step_rel;