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"
72 #include "domdec_network.h"
78 #include "compute_io.h"
80 #include "checkpoint.h"
81 #include "mtop_util.h"
82 #include "sighandler.h"
85 #include "pme_loadbal.h"
88 #include "types/nlistheuristics.h"
89 #include "types/iteratedconstraints.h"
90 #include "nbnxn_cuda_data_mgmt.h"
92 #include "gromacs/utility/gmxmpi.h"
98 static void reset_all_counters(FILE *fplog, t_commrec *cr,
100 gmx_large_int_t *step_rel, t_inputrec *ir,
101 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
102 gmx_runtime_t *runtime,
103 nbnxn_cuda_ptr_t cu_nbv)
105 char sbuf[STEPSTRSIZE];
107 /* Reset all the counters related to performance over the run */
108 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
109 gmx_step_str(step, sbuf));
113 nbnxn_cuda_reset_timings(cu_nbv);
116 wallcycle_stop(wcycle, ewcRUN);
117 wallcycle_reset_all(wcycle);
118 if (DOMAINDECOMP(cr))
120 reset_dd_statistics_counters(cr->dd);
123 ir->init_step += *step_rel;
124 ir->nsteps -= *step_rel;
126 wallcycle_start(wcycle, ewcRUN);
127 runtime_start(runtime);
128 print_date_and_time(fplog, cr->nodeid, "Restarted time", runtime);
131 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
132 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
134 gmx_vsite_t *vsite, gmx_constr_t constr,
135 int stepout, t_inputrec *ir,
136 gmx_mtop_t *top_global,
138 t_state *state_global,
140 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
141 gmx_edsam_t ed, t_forcerec *fr,
142 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
143 real cpt_period, real max_hours,
144 const char gmx_unused *deviceOptions,
146 gmx_runtime_t *runtime)
149 gmx_large_int_t step, step_rel;
151 double t, t0, lam0[efptNR];
152 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
153 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
154 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
155 bBornRadii, bStartingFromCpt;
156 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
157 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
158 bForceUpdate = FALSE, bCPT;
160 gmx_bool bMasterState;
161 int force_flags, cglo_flags;
162 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
167 t_state *bufstate = NULL;
168 matrix *scale_tot, pcoupl_mu, M, ebox;
171 gmx_repl_ex_t repl_ex = NULL;
174 t_mdebin *mdebin = NULL;
175 df_history_t df_history;
176 t_state *state = NULL;
177 rvec *f_global = NULL;
180 gmx_enerdata_t *enerd;
182 gmx_global_stat_t gstat;
183 gmx_update_t upd = NULL;
184 t_graph *graph = NULL;
186 gmx_rng_t mcrng = NULL;
187 gmx_groups_t *groups;
188 gmx_ekindata_t *ekind, *ekind_save;
189 gmx_shellfc_t shellfc;
190 int count, nconverged = 0;
193 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
195 gmx_bool bResetCountersHalfMaxH = FALSE;
196 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
197 gmx_bool bUpdateDoLR;
202 real veta_save, scalevir, tracevir;
208 real saved_conserved_quantity = 0;
213 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
214 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
215 gmx_iterate_t iterate;
216 gmx_large_int_t multisim_nsteps = -1; /* number of steps to do before first multisim
217 simulation stops. If equal to zero, don't
218 communicate any more between multisims.*/
219 /* PME load balancing data for GPU kernels */
220 pme_load_balancing_t pme_loadbal = NULL;
222 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
225 /* Temporary addition for FAHCORE checkpointing */
229 /* Check for special mdrun options */
230 bRerunMD = (Flags & MD_RERUN);
231 bAppend = (Flags & MD_APPENDFILES);
232 if (Flags & MD_RESETCOUNTERSHALFWAY)
236 /* Signal to reset the counters half the simulation steps. */
237 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
239 /* Signal to reset the counters halfway the simulation time. */
240 bResetCountersHalfMaxH = (max_hours > 0);
243 /* md-vv uses averaged full step velocities for T-control
244 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
245 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
247 if (bVV) /* to store the initial velocities while computing virial */
249 snew(cbuf, top_global->natoms);
251 /* all the iteratative cases - only if there are constraints */
252 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
253 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
254 false in this step. The correct value, true or false,
255 is set at each step, as it depends on the frequency of temperature
256 and pressure control.*/
257 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
261 /* Since we don't know if the frames read are related in any way,
262 * rebuild the neighborlist at every step.
265 ir->nstcalcenergy = 1;
269 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
271 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
272 bGStatEveryStep = (nstglobalcomm == 1);
274 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
277 "To reduce the energy communication with nstlist = -1\n"
278 "the neighbor list validity should not be checked at every step,\n"
279 "this means that exact integration is not guaranteed.\n"
280 "The neighbor list validity is checked after:\n"
281 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
282 "In most cases this will result in exact integration.\n"
283 "This reduces the energy communication by a factor of 2 to 3.\n"
284 "If you want less energy communication, set nstlist > 3.\n\n");
291 groups = &top_global->groups;
294 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
295 &(state_global->fep_state), lam0,
296 nrnb, top_global, &upd,
297 nfile, fnm, &outf, &mdebin,
298 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
300 clear_mat(total_vir);
302 /* Energy terms and groups */
304 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
306 if (DOMAINDECOMP(cr))
312 snew(f, top_global->natoms);
315 /* lambda Monte carlo random number generator */
318 mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
320 /* copy the state into df_history */
321 copy_df_history(&df_history, &state_global->dfhist);
323 /* Kinetic energy data */
325 init_ekindata(fplog, top_global, &(ir->opts), ekind);
326 /* needed for iteration of constraints */
328 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
329 /* Copy the cos acceleration to the groups struct */
330 ekind->cosacc.cos_accel = ir->cos_accel;
332 gstat = global_stat_init(ir);
335 /* Check for polarizable models and flexible constraints */
336 shellfc = init_shell_flexcon(fplog,
337 top_global, n_flexible_constraints(constr),
338 (ir->bContinuation ||
339 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
340 NULL : state_global->x);
344 #ifdef GMX_THREAD_MPI
345 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
347 set_deform_reference_box(upd,
348 deform_init_init_step_tpx,
349 deform_init_box_tpx);
350 #ifdef GMX_THREAD_MPI
351 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
356 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
357 if ((io > 2000) && MASTER(cr))
360 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
365 if (DOMAINDECOMP(cr))
367 top = dd_init_local_top(top_global);
370 dd_init_local_state(cr->dd, state_global, state);
372 if (DDMASTER(cr->dd) && ir->nstfout)
374 snew(f_global, state_global->natoms);
381 /* Initialize the particle decomposition and split the topology */
382 top = split_system(fplog, top_global, ir, cr);
384 pd_cg_range(cr, &fr->cg0, &fr->hcg);
385 pd_at_range(cr, &a0, &a1);
389 top = gmx_mtop_generate_local_top(top_global, ir);
392 a1 = top_global->natoms;
395 forcerec_set_excl_load(fr, top, cr);
397 state = partdec_init_local_state(cr, state_global);
400 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
404 set_vsite_top(vsite, top, mdatoms, cr);
407 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
409 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
414 make_local_shells(cr, mdatoms, shellfc);
417 setup_bonded_threading(fr, &top->idef);
419 if (ir->pull && PAR(cr))
421 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
425 if (DOMAINDECOMP(cr))
427 /* Distribute the charge groups over the nodes from the master node */
428 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
429 state_global, top_global, ir,
430 state, &f, mdatoms, top, fr,
431 vsite, shellfc, constr,
432 nrnb, wcycle, FALSE);
436 update_mdatoms(mdatoms, state->lambda[efptMASS]);
438 if (opt2bSet("-cpi", nfile, fnm))
440 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
444 bStateFromCP = FALSE;
451 /* Update mdebin with energy history if appending to output files */
452 if (Flags & MD_APPENDFILES)
454 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
458 /* We might have read an energy history from checkpoint,
459 * free the allocated memory and reset the counts.
461 done_energyhistory(&state_global->enerhist);
462 init_energyhistory(&state_global->enerhist);
465 /* Set the initial energy history in state by updating once */
466 update_energyhistory(&state_global->enerhist, mdebin);
469 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
471 /* Set the random state if we read a checkpoint file */
472 set_stochd_state(upd, state);
475 if (state->flags & (1<<estMC_RNG))
477 set_mc_state(mcrng, state);
480 /* Initialize constraints */
483 if (!DOMAINDECOMP(cr))
485 set_constraints(constr, top, ir, mdatoms, cr);
491 /* We need to be sure replica exchange can only occur
492 * when the energies are current */
493 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
494 "repl_ex_nst", &repl_ex_nst);
495 /* This check needs to happen before inter-simulation
496 * signals are initialized, too */
498 if (repl_ex_nst > 0 && MASTER(cr))
500 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
501 repl_ex_nst, repl_ex_nex, repl_ex_seed);
504 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
505 * With perturbed charges with soft-core we should not change the cut-off.
507 if ((Flags & MD_TUNEPME) &&
508 EEL_PME(fr->eeltype) &&
509 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
510 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
513 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
515 if (cr->duty & DUTY_PME)
517 /* Start tuning right away, as we can't measure the load */
518 bPMETuneRunning = TRUE;
522 /* Separate PME nodes, we can measure the PP/PME load balance */
527 if (!ir->bContinuation && !bRerunMD)
529 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
531 /* Set the velocities of frozen particles to zero */
532 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
534 for (m = 0; m < DIM; m++)
536 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
546 /* Constrain the initial coordinates and velocities */
547 do_constrain_first(fplog, constr, ir, mdatoms, state,
552 /* Construct the virtual sites for the initial configuration */
553 construct_vsites(vsite, state->x, ir->delta_t, NULL,
554 top->idef.iparams, top->idef.il,
555 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
561 /* set free energy calculation frequency as the minimum
562 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
563 nstfep = ir->fepvals->nstdhdl;
566 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
570 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
573 /* I'm assuming we need global communication the first time! MRS */
574 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
575 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
576 | (bVV ? CGLO_PRESSURE : 0)
577 | (bVV ? CGLO_CONSTRAINT : 0)
578 | (bRerunMD ? CGLO_RERUNMD : 0)
579 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
581 bSumEkinhOld = FALSE;
582 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
583 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
584 constr, NULL, FALSE, state->box,
585 top_global, &bSumEkinhOld, cglo_flags);
586 if (ir->eI == eiVVAK)
588 /* a second call to get the half step temperature initialized as well */
589 /* we do the same call as above, but turn the pressure off -- internally to
590 compute_globals, this is recognized as a velocity verlet half-step
591 kinetic energy calculation. This minimized excess variables, but
592 perhaps loses some logic?*/
594 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
595 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
596 constr, NULL, FALSE, state->box,
597 top_global, &bSumEkinhOld,
598 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
601 /* Calculate the initial half step temperature, and save the ekinh_old */
602 if (!(Flags & MD_STARTFROMCPT))
604 for (i = 0; (i < ir->opts.ngtc); i++)
606 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
611 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
612 and there is no previous step */
615 /* if using an iterative algorithm, we need to create a working directory for the state. */
618 bufstate = init_bufstate(state);
621 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
622 temperature control */
623 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
627 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
630 "RMS relative constraint deviation after constraining: %.2e\n",
631 constr_rmsd(constr, FALSE));
633 if (EI_STATE_VELOCITY(ir->eI))
635 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
639 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
640 " input trajectory '%s'\n\n",
641 *(top_global->name), opt2fn("-rerun", nfile, fnm));
644 fprintf(stderr, "Calculated time to finish depends on nsteps from "
645 "run input file,\nwhich may not correspond to the time "
646 "needed to process input trajectory.\n\n");
652 fprintf(stderr, "starting mdrun '%s'\n",
653 *(top_global->name));
656 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
660 sprintf(tbuf, "%s", "infinite");
662 if (ir->init_step > 0)
664 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
665 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
666 gmx_step_str(ir->init_step, sbuf2),
667 ir->init_step*ir->delta_t);
671 fprintf(stderr, "%s steps, %s ps.\n",
672 gmx_step_str(ir->nsteps, sbuf), tbuf);
675 fprintf(fplog, "\n");
678 /* Set and write start time */
679 runtime_start(runtime);
680 print_date_and_time(fplog, cr->nodeid, "Started mdrun", runtime);
681 wallcycle_start(wcycle, ewcRUN);
684 fprintf(fplog, "\n");
687 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
689 chkpt_ret = fcCheckPointParallel( cr->nodeid,
693 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
698 /***********************************************************
702 ************************************************************/
704 /* if rerunMD then read coordinates and velocities from input trajectory */
707 if (getenv("GMX_FORCE_UPDATE"))
715 bNotLastFrame = read_first_frame(oenv, &status,
716 opt2fn("-rerun", nfile, fnm),
717 &rerun_fr, TRX_NEED_X | TRX_READ_V);
718 if (rerun_fr.natoms != top_global->natoms)
721 "Number of atoms in trajectory (%d) does not match the "
722 "run input file (%d)\n",
723 rerun_fr.natoms, top_global->natoms);
725 if (ir->ePBC != epbcNONE)
729 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);
731 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
733 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
740 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
743 if (ir->ePBC != epbcNONE)
745 /* Set the shift vectors.
746 * Necessary here when have a static box different from the tpr box.
748 calc_shifts(rerun_fr.box, fr->shift_vec);
752 /* loop over MD steps or if rerunMD to end of input trajectory */
754 /* Skip the first Nose-Hoover integration when we get the state from tpx */
755 bStateFromTPX = !bStateFromCP;
756 bInitStep = bFirstStep && (bStateFromTPX || bVV);
757 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
759 bSumEkinhOld = FALSE;
762 init_global_signals(&gs, cr, ir, repl_ex_nst);
764 step = ir->init_step;
767 if (ir->nstlist == -1)
769 init_nlistheuristics(&nlh, bGStatEveryStep, step);
772 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
774 /* check how many steps are left in other sims */
775 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
779 /* and stop now if we should */
780 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
781 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
782 while (!bLastStep || (bRerunMD && bNotLastFrame))
785 wallcycle_start(wcycle, ewcSTEP);
791 step = rerun_fr.step;
792 step_rel = step - ir->init_step;
805 bLastStep = (step_rel == ir->nsteps);
806 t = t0 + step*ir->delta_t;
809 if (ir->efep != efepNO || ir->bSimTemp)
811 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
812 requiring different logic. */
814 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
815 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
816 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
817 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded) && (step > 0));
822 update_annealing_target_temp(&(ir->opts), t);
827 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
829 for (i = 0; i < state_global->natoms; i++)
831 copy_rvec(rerun_fr.x[i], state_global->x[i]);
835 for (i = 0; i < state_global->natoms; i++)
837 copy_rvec(rerun_fr.v[i], state_global->v[i]);
842 for (i = 0; i < state_global->natoms; i++)
844 clear_rvec(state_global->v[i]);
848 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
849 " Ekin, temperature and pressure are incorrect,\n"
850 " the virial will be incorrect when constraints are present.\n"
852 bRerunWarnNoV = FALSE;
856 copy_mat(rerun_fr.box, state_global->box);
857 copy_mat(state_global->box, state->box);
859 if (vsite && (Flags & MD_RERUN_VSITE))
861 if (DOMAINDECOMP(cr))
863 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
867 /* Following is necessary because the graph may get out of sync
868 * with the coordinates if we only have every N'th coordinate set
870 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
871 shift_self(graph, state->box, state->x);
873 construct_vsites(vsite, state->x, ir->delta_t, state->v,
874 top->idef.iparams, top->idef.il,
875 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
878 unshift_self(graph, state->box, state->x);
883 /* Stop Center of Mass motion */
884 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
888 /* for rerun MD always do Neighbour Searching */
889 bNS = (bFirstStep || ir->nstlist != 0);
894 /* Determine whether or not to do Neighbour Searching and LR */
895 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
897 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
898 (ir->nstlist == -1 && nlh.nabnsb > 0));
900 if (bNS && ir->nstlist == -1)
902 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
906 /* check whether we should stop because another simulation has
910 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
911 (multisim_nsteps != ir->nsteps) )
918 "Stopping simulation %d because another one has finished\n",
922 gs.sig[eglsCHKPT] = 1;
927 /* < 0 means stop at next step, > 0 means stop at next NS step */
928 if ( (gs.set[eglsSTOPCOND] < 0) ||
929 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
934 /* Determine whether or not to update the Born radii if doing GB */
935 bBornRadii = bFirstStep;
936 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
941 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
942 do_verbose = bVerbose &&
943 (step % stepout == 0 || bFirstStep || bLastStep);
945 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
953 bMasterState = FALSE;
954 /* Correct the new box if it is too skewed */
955 if (DYNAMIC_BOX(*ir))
957 if (correct_box(fplog, step, state->box, graph))
962 if (DOMAINDECOMP(cr) && bMasterState)
964 dd_collect_state(cr->dd, state, state_global);
968 if (DOMAINDECOMP(cr))
970 /* Repartition the domain decomposition */
971 wallcycle_start(wcycle, ewcDOMDEC);
972 dd_partition_system(fplog, step, cr,
973 bMasterState, nstglobalcomm,
974 state_global, top_global, ir,
975 state, &f, mdatoms, top, fr,
976 vsite, shellfc, constr,
978 do_verbose && !bPMETuneRunning);
979 wallcycle_stop(wcycle, ewcDOMDEC);
980 /* If using an iterative integrator, reallocate space to match the decomposition */
984 if (MASTER(cr) && do_log)
986 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
989 if (ir->efep != efepNO)
991 update_mdatoms(mdatoms, state->lambda[efptMASS]);
994 if ((bRerunMD && rerun_fr.bV) || bExchanged)
997 /* We need the kinetic energy at minus the half step for determining
998 * the full step kinetic energy and possibly for T-coupling.*/
999 /* This may not be quite working correctly yet . . . . */
1000 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1001 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1002 constr, NULL, FALSE, state->box,
1003 top_global, &bSumEkinhOld,
1004 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1006 clear_mat(force_vir);
1008 /* We write a checkpoint at this MD step when:
1009 * either at an NS step when we signalled through gs,
1010 * or at the last step (but not when we do not want confout),
1011 * but never at the first step or with rerun.
1013 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1014 (bLastStep && (Flags & MD_CONFOUT))) &&
1015 step > ir->init_step && !bRerunMD);
1018 gs.set[eglsCHKPT] = 0;
1021 /* Determine the energy and pressure:
1022 * at nstcalcenergy steps and at energy output steps (set below).
1024 if (EI_VV(ir->eI) && (!bInitStep))
1026 /* for vv, the first half of the integration actually corresponds
1027 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1028 but the virial needs to be calculated on both the current step and the 'next' step. Future
1029 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1031 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1032 bCalcVir = bCalcEner ||
1033 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1037 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1038 bCalcVir = bCalcEner ||
1039 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1042 /* Do we need global communication ? */
1043 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1044 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1045 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1047 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1049 if (do_ene || do_log)
1056 /* these CGLO_ options remain the same throughout the iteration */
1057 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1058 (bGStat ? CGLO_GSTAT : 0)
1061 force_flags = (GMX_FORCE_STATECHANGED |
1062 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1063 GMX_FORCE_ALLFORCES |
1065 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1066 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1067 (bDoFEP ? GMX_FORCE_DHDL : 0)
1072 if (do_per_step(step, ir->nstcalclr))
1074 force_flags |= GMX_FORCE_DO_LR;
1080 /* Now is the time to relax the shells */
1081 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1082 ir, bNS, force_flags,
1085 state, f, force_vir, mdatoms,
1086 nrnb, wcycle, graph, groups,
1087 shellfc, fr, bBornRadii, t, mu_tot,
1099 /* The coordinates (x) are shifted (to get whole molecules)
1101 * This is parallellized as well, and does communication too.
1102 * Check comments in sim_util.c
1104 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1105 state->box, state->x, &state->hist,
1106 f, force_vir, mdatoms, enerd, fcd,
1107 state->lambda, graph,
1108 fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1109 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1112 if (bVV && !bStartingFromCpt && !bRerunMD)
1113 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1115 if (ir->eI == eiVV && bInitStep)
1117 /* if using velocity verlet with full time step Ekin,
1118 * take the first half step only to compute the
1119 * virial for the first step. From there,
1120 * revert back to the initial coordinates
1121 * so that the input is actually the initial step.
1123 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1127 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1128 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1131 /* If we are using twin-range interactions where the long-range component
1132 * is only evaluated every nstcalclr>1 steps, we should do a special update
1133 * step to combine the long-range forces on these steps.
1134 * For nstcalclr=1 this is not done, since the forces would have been added
1135 * directly to the short-range forces already.
1137 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1139 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1140 f, bUpdateDoLR, fr->f_twin, fcd,
1141 ekind, M, upd, bInitStep, etrtVELOCITY1,
1142 cr, nrnb, constr, &top->idef);
1144 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1146 gmx_iterate_init(&iterate, TRUE);
1148 /* for iterations, we save these vectors, as we will be self-consistently iterating
1151 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1153 /* save the state */
1154 if (iterate.bIterationActive)
1156 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1159 bFirstIterate = TRUE;
1160 while (bFirstIterate || iterate.bIterationActive)
1162 if (iterate.bIterationActive)
1164 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1165 if (bFirstIterate && bTrotter)
1167 /* The first time through, we need a decent first estimate
1168 of veta(t+dt) to compute the constraints. Do
1169 this by computing the box volume part of the
1170 trotter integration at this time. Nothing else
1171 should be changed by this routine here. If
1172 !(first time), we start with the previous value
1175 veta_save = state->veta;
1176 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1177 vetanew = state->veta;
1178 state->veta = veta_save;
1183 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1185 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1186 state, fr->bMolPBC, graph, f,
1187 &top->idef, shake_vir,
1188 cr, nrnb, wcycle, upd, constr,
1189 TRUE, bCalcVir, vetanew);
1193 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1199 /* Need to unshift here if a do_force has been
1200 called in the previous step */
1201 unshift_self(graph, state->box, state->x);
1204 /* if VV, compute the pressure and constraints */
1205 /* For VV2, we strictly only need this if using pressure
1206 * control, but we really would like to have accurate pressures
1208 * Think about ways around this in the future?
1209 * For now, keep this choice in comments.
1211 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1212 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1214 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1215 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1217 bSumEkinhOld = TRUE;
1219 /* for vv, the first half of the integration actually corresponds to the previous step.
1220 So we need information from the last step in the first half of the integration */
1221 if (bGStat || do_per_step(step-1, nstglobalcomm))
1223 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1224 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1225 constr, NULL, FALSE, state->box,
1226 top_global, &bSumEkinhOld,
1229 | (bTemp ? CGLO_TEMPERATURE : 0)
1230 | (bPres ? CGLO_PRESSURE : 0)
1231 | (bPres ? CGLO_CONSTRAINT : 0)
1232 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1233 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1236 /* explanation of above:
1237 a) We compute Ekin at the full time step
1238 if 1) we are using the AveVel Ekin, and it's not the
1239 initial step, or 2) if we are using AveEkin, but need the full
1240 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1241 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1242 EkinAveVel because it's needed for the pressure */
1244 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1249 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1250 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1257 /* We need the kinetic energy at minus the half step for determining
1258 * the full step kinetic energy and possibly for T-coupling.*/
1259 /* This may not be quite working correctly yet . . . . */
1260 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1261 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1262 constr, NULL, FALSE, state->box,
1263 top_global, &bSumEkinhOld,
1264 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1269 if (iterate.bIterationActive &&
1270 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1271 state->veta, &vetanew))
1275 bFirstIterate = FALSE;
1278 if (bTrotter && !bInitStep)
1280 copy_mat(shake_vir, state->svir_prev);
1281 copy_mat(force_vir, state->fvir_prev);
1282 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1284 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1285 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1286 enerd->term[F_EKIN] = trace(ekind->ekin);
1289 /* if it's the initial step, we performed this first step just to get the constraint virial */
1290 if (bInitStep && ir->eI == eiVV)
1292 copy_rvecn(cbuf, state->v, 0, state->natoms);
1296 /* MRS -- now done iterating -- compute the conserved quantity */
1299 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1302 last_ekin = enerd->term[F_EKIN];
1304 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1306 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1308 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1311 sum_dhdl(enerd, state->lambda, ir->fepvals);
1315 /* ######## END FIRST UPDATE STEP ############## */
1316 /* ######## If doing VV, we now have v(dt) ###### */
1319 /* perform extended ensemble sampling in lambda - we don't
1320 actually move to the new state before outputting
1321 statistics, but if performing simulated tempering, we
1322 do update the velocities and the tau_t. */
1324 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, &df_history, step, mcrng, state->v, mdatoms);
1326 /* ################## START TRAJECTORY OUTPUT ################# */
1328 /* Now we have the energies and forces corresponding to the
1329 * coordinates at time t. We must output all of this before
1331 * for RerunMD t is read from input trajectory
1334 if (do_per_step(step, ir->nstxout))
1336 mdof_flags |= MDOF_X;
1338 if (do_per_step(step, ir->nstvout))
1340 mdof_flags |= MDOF_V;
1342 if (do_per_step(step, ir->nstfout))
1344 mdof_flags |= MDOF_F;
1346 if (do_per_step(step, ir->nstxtcout))
1348 mdof_flags |= MDOF_XTC;
1352 mdof_flags |= MDOF_CPT;
1356 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1359 /* Enforce writing positions and velocities at end of run */
1360 mdof_flags |= (MDOF_X | MDOF_V);
1366 fcReportProgress( ir->nsteps, step );
1369 /* sync bCPT and fc record-keeping */
1370 if (bCPT && MASTER(cr))
1372 fcRequestCheckPoint();
1376 if (mdof_flags != 0)
1378 wallcycle_start(wcycle, ewcTRAJ);
1381 if (state->flags & (1<<estLD_RNG))
1383 get_stochd_state(upd, state);
1385 if (state->flags & (1<<estMC_RNG))
1387 get_mc_state(mcrng, state);
1393 state_global->ekinstate.bUpToDate = FALSE;
1397 update_ekinstate(&state_global->ekinstate, ekind);
1398 state_global->ekinstate.bUpToDate = TRUE;
1400 update_energyhistory(&state_global->enerhist, mdebin);
1401 if (ir->efep != efepNO || ir->bSimTemp)
1403 state_global->fep_state = state->fep_state; /* MRS: seems kludgy. The code should be
1404 structured so this isn't necessary.
1405 Note this reassignment is only necessary
1406 for single threads.*/
1407 copy_df_history(&state_global->dfhist, &df_history);
1411 write_traj(fplog, cr, outf, mdof_flags, top_global,
1412 step, t, state, state_global, f, f_global, &n_xtc, &x_xtc);
1419 if (bLastStep && step_rel == ir->nsteps &&
1420 (Flags & MD_CONFOUT) && MASTER(cr) &&
1423 /* x and v have been collected in write_traj,
1424 * because a checkpoint file will always be written
1427 fprintf(stderr, "\nWriting final coordinates.\n");
1430 /* Make molecules whole only for confout writing */
1431 do_pbc_mtop(fplog, ir->ePBC, state->box, top_global, state_global->x);
1433 write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm),
1434 *top_global->name, top_global,
1435 state_global->x, state_global->v,
1436 ir->ePBC, state->box);
1439 wallcycle_stop(wcycle, ewcTRAJ);
1442 /* kludge -- virial is lost with restart for NPT control. Must restart */
1443 if (bStartingFromCpt && bVV)
1445 copy_mat(state->svir_prev, shake_vir);
1446 copy_mat(state->fvir_prev, force_vir);
1448 /* ################## END TRAJECTORY OUTPUT ################ */
1450 /* Determine the wallclock run time up till now */
1451 run_time = gmx_gettime() - (double)runtime->real;
1453 /* Check whether everything is still allright */
1454 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1455 #ifdef GMX_THREAD_MPI
1460 /* this is just make gs.sig compatible with the hack
1461 of sending signals around by MPI_Reduce with together with
1463 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1465 gs.sig[eglsSTOPCOND] = 1;
1467 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1469 gs.sig[eglsSTOPCOND] = -1;
1471 /* < 0 means stop at next step, > 0 means stop at next NS step */
1475 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1476 gmx_get_signal_name(),
1477 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1481 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1482 gmx_get_signal_name(),
1483 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1485 handled_stop_condition = (int)gmx_get_stop_condition();
1487 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1488 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1489 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1491 /* Signal to terminate the run */
1492 gs.sig[eglsSTOPCOND] = 1;
1495 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1497 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1500 if (bResetCountersHalfMaxH && MASTER(cr) &&
1501 run_time > max_hours*60.0*60.0*0.495)
1503 gs.sig[eglsRESETCOUNTERS] = 1;
1506 if (ir->nstlist == -1 && !bRerunMD)
1508 /* When bGStatEveryStep=FALSE, global_stat is only called
1509 * when we check the atom displacements, not at NS steps.
1510 * This means that also the bonded interaction count check is not
1511 * performed immediately after NS. Therefore a few MD steps could
1512 * be performed with missing interactions.
1513 * But wrong energies are never written to file,
1514 * since energies are only written after global_stat
1517 if (step >= nlh.step_nscheck)
1519 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1520 nlh.scale_tot, state->x);
1524 /* This is not necessarily true,
1525 * but step_nscheck is determined quite conservatively.
1531 /* In parallel we only have to check for checkpointing in steps
1532 * where we do global communication,
1533 * otherwise the other nodes don't know.
1535 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1538 run_time >= nchkpt*cpt_period*60.0)) &&
1539 gs.set[eglsCHKPT] == 0)
1541 gs.sig[eglsCHKPT] = 1;
1544 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1549 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1551 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1553 gmx_bool bIfRandomize;
1554 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1555 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1556 if (constr && bIfRandomize)
1558 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1559 state, fr->bMolPBC, graph, f,
1560 &top->idef, tmp_vir,
1561 cr, nrnb, wcycle, upd, constr,
1562 TRUE, bCalcVir, vetanew);
1567 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1569 gmx_iterate_init(&iterate, TRUE);
1570 /* for iterations, we save these vectors, as we will be redoing the calculations */
1571 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1574 bFirstIterate = TRUE;
1575 while (bFirstIterate || iterate.bIterationActive)
1577 /* We now restore these vectors to redo the calculation with improved extended variables */
1578 if (iterate.bIterationActive)
1580 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1583 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1584 so scroll down for that logic */
1586 /* ######### START SECOND UPDATE STEP ################# */
1587 /* Box is changed in update() when we do pressure coupling,
1588 * but we should still use the old box for energy corrections and when
1589 * writing it to the energy file, so it matches the trajectory files for
1590 * the same timestep above. Make a copy in a separate array.
1592 copy_mat(state->box, lastbox);
1597 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1599 wallcycle_start(wcycle, ewcUPDATE);
1600 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1603 if (iterate.bIterationActive)
1611 /* we use a new value of scalevir to converge the iterations faster */
1612 scalevir = tracevir/trace(shake_vir);
1614 msmul(shake_vir, scalevir, shake_vir);
1615 m_add(force_vir, shake_vir, total_vir);
1616 clear_mat(shake_vir);
1618 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1619 /* We can only do Berendsen coupling after we have summed
1620 * the kinetic energy or virial. Since the happens
1621 * in global_state after update, we should only do it at
1622 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1627 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1628 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1633 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1635 /* velocity half-step update */
1636 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1637 bUpdateDoLR, fr->f_twin, fcd,
1638 ekind, M, upd, FALSE, etrtVELOCITY2,
1639 cr, nrnb, constr, &top->idef);
1642 /* Above, initialize just copies ekinh into ekin,
1643 * it doesn't copy position (for VV),
1644 * and entire integrator for MD.
1647 if (ir->eI == eiVVAK)
1649 copy_rvecn(state->x, cbuf, 0, state->natoms);
1651 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1653 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1654 bUpdateDoLR, fr->f_twin, fcd,
1655 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1656 wallcycle_stop(wcycle, ewcUPDATE);
1658 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1659 fr->bMolPBC, graph, f,
1660 &top->idef, shake_vir,
1661 cr, nrnb, wcycle, upd, constr,
1662 FALSE, bCalcVir, state->veta);
1664 if (ir->eI == eiVVAK)
1666 /* erase F_EKIN and F_TEMP here? */
1667 /* just compute the kinetic energy at the half step to perform a trotter step */
1668 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1669 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1670 constr, NULL, FALSE, lastbox,
1671 top_global, &bSumEkinhOld,
1672 cglo_flags | CGLO_TEMPERATURE
1674 wallcycle_start(wcycle, ewcUPDATE);
1675 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1676 /* now we know the scaling, we can compute the positions again again */
1677 copy_rvecn(cbuf, state->x, 0, state->natoms);
1679 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1681 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1682 bUpdateDoLR, fr->f_twin, fcd,
1683 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1684 wallcycle_stop(wcycle, ewcUPDATE);
1686 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1687 /* are the small terms in the shake_vir here due
1688 * to numerical errors, or are they important
1689 * physically? I'm thinking they are just errors, but not completely sure.
1690 * For now, will call without actually constraining, constr=NULL*/
1691 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1692 state, fr->bMolPBC, graph, f,
1693 &top->idef, tmp_vir,
1694 cr, nrnb, wcycle, upd, NULL,
1700 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1703 if (fr->bSepDVDL && fplog && do_log)
1705 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1709 /* this factor or 2 correction is necessary
1710 because half of the constraint force is removed
1711 in the vv step, so we have to double it. See
1712 the Redmine issue #1255. It is not yet clear
1713 if the factor of 2 is exact, or just a very
1714 good approximation, and this will be
1715 investigated. The next step is to see if this
1716 can be done adding a dhdl contribution from the
1717 rattle step, but this is somewhat more
1718 complicated with the current code. Will be
1719 investigated, hopefully for 4.6.3. However,
1720 this current solution is much better than
1721 having it completely wrong.
1723 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1727 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1732 /* Need to unshift here */
1733 unshift_self(graph, state->box, state->x);
1738 wallcycle_start(wcycle, ewcVSITECONSTR);
1741 shift_self(graph, state->box, state->x);
1743 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1744 top->idef.iparams, top->idef.il,
1745 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1749 unshift_self(graph, state->box, state->x);
1751 wallcycle_stop(wcycle, ewcVSITECONSTR);
1754 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1755 /* With Leap-Frog we can skip compute_globals at
1756 * non-communication steps, but we need to calculate
1757 * the kinetic energy one step before communication.
1759 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1761 if (ir->nstlist == -1 && bFirstIterate)
1763 gs.sig[eglsNABNSB] = nlh.nabnsb;
1765 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1766 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1768 bFirstIterate ? &gs : NULL,
1769 (step_rel % gs.nstms == 0) &&
1770 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1772 top_global, &bSumEkinhOld,
1774 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1775 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1776 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1777 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1778 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1779 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1782 if (ir->nstlist == -1 && bFirstIterate)
1784 nlh.nabnsb = gs.set[eglsNABNSB];
1785 gs.set[eglsNABNSB] = 0;
1788 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1789 /* ############# END CALC EKIN AND PRESSURE ################# */
1791 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1792 the virial that should probably be addressed eventually. state->veta has better properies,
1793 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1794 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1796 if (iterate.bIterationActive &&
1797 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1798 trace(shake_vir), &tracevir))
1802 bFirstIterate = FALSE;
1805 if (!bVV || bRerunMD)
1807 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1808 sum_dhdl(enerd, state->lambda, ir->fepvals);
1810 update_box(fplog, step, ir, mdatoms, state, f,
1811 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1813 /* ################# END UPDATE STEP 2 ################# */
1814 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1816 /* The coordinates (x) were unshifted in update */
1819 /* We will not sum ekinh_old,
1820 * so signal that we still have to do it.
1822 bSumEkinhOld = TRUE;
1825 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1827 /* use the directly determined last velocity, not actually the averaged half steps */
1828 if (bTrotter && ir->eI == eiVV)
1830 enerd->term[F_EKIN] = last_ekin;
1832 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1836 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1840 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1842 /* ######### END PREPARING EDR OUTPUT ########### */
1844 /* Time for performance */
1845 if (((step % stepout) == 0) || bLastStep)
1847 runtime_upd_proc(runtime);
1853 gmx_bool do_dr, do_or;
1855 if (fplog && do_log && bDoExpanded)
1857 /* only needed if doing expanded ensemble */
1858 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1859 &df_history, state->fep_state, ir->nstlog, step);
1861 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1865 upd_mdebin(mdebin, bDoDHDL, TRUE,
1866 t, mdatoms->tmass, enerd, state,
1867 ir->fepvals, ir->expandedvals, lastbox,
1868 shake_vir, force_vir, total_vir, pres,
1869 ekind, mu_tot, constr);
1873 upd_mdebin_step(mdebin);
1876 do_dr = do_per_step(step, ir->nstdisreout);
1877 do_or = do_per_step(step, ir->nstorireout);
1879 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
1881 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1883 if (ir->ePull != epullNO)
1885 pull_print_output(ir->pull, step, t);
1888 if (do_per_step(step, ir->nstlog))
1890 if (fflush(fplog) != 0)
1892 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1898 /* Have to do this part after outputting the logfile and the edr file */
1899 state->fep_state = lamnew;
1900 for (i = 0; i < efptNR; i++)
1902 state_global->lambda[i] = ir->fepvals->all_lambda[i][lamnew];
1905 /* Remaining runtime */
1906 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1910 fprintf(stderr, "\n");
1912 print_time(stderr, runtime, step, ir, cr);
1915 /* Replica exchange */
1917 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1918 do_per_step(step, repl_ex_nst))
1920 bExchanged = replica_exchange(fplog, cr, repl_ex,
1921 state_global, enerd,
1924 if (bExchanged && DOMAINDECOMP(cr))
1926 dd_partition_system(fplog, step, cr, TRUE, 1,
1927 state_global, top_global, ir,
1928 state, &f, mdatoms, top, fr,
1929 vsite, shellfc, constr,
1930 nrnb, wcycle, FALSE);
1936 bStartingFromCpt = FALSE;
1938 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1939 /* With all integrators, except VV, we need to retain the pressure
1940 * at the current step for coupling at the next step.
1942 if ((state->flags & (1<<estPRES_PREV)) &&
1944 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1946 /* Store the pressure in t_state for pressure coupling
1947 * at the next MD step.
1949 copy_mat(pres, state->pres_prev);
1952 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1954 if ( (membed != NULL) && (!bLastStep) )
1956 rescale_membed(step_rel, membed, state_global->x);
1963 /* read next frame from input trajectory */
1964 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1969 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1973 if (!bRerunMD || !rerun_fr.bStep)
1975 /* increase the MD step number */
1980 cycles = wallcycle_stop(wcycle, ewcSTEP);
1981 if (DOMAINDECOMP(cr) && wcycle)
1983 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1986 if (bPMETuneRunning || bPMETuneTry)
1988 /* PME grid + cut-off optimization with GPUs or PME nodes */
1990 /* Count the total cycles over the last steps */
1991 cycles_pmes += cycles;
1993 /* We can only switch cut-off at NS steps */
1994 if (step % ir->nstlist == 0)
1996 /* PME grid + cut-off optimization with GPUs or PME nodes */
1999 if (DDMASTER(cr->dd))
2001 /* PME node load is too high, start tuning */
2002 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
2004 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
2006 if (bPMETuneRunning || step_rel > ir->nstlist*50)
2008 bPMETuneTry = FALSE;
2011 if (bPMETuneRunning)
2013 /* init_step might not be a multiple of nstlist,
2014 * but the first cycle is always skipped anyhow.
2017 pme_load_balance(pme_loadbal, cr,
2018 (bVerbose && MASTER(cr)) ? stderr : NULL,
2020 ir, state, cycles_pmes,
2021 fr->ic, fr->nbv, &fr->pmedata,
2024 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2025 fr->ewaldcoeff = fr->ic->ewaldcoeff;
2026 fr->rlist = fr->ic->rlist;
2027 fr->rlistlong = fr->ic->rlistlong;
2028 fr->rcoulomb = fr->ic->rcoulomb;
2029 fr->rvdw = fr->ic->rvdw;
2035 if (step_rel == wcycle_get_reset_counters(wcycle) ||
2036 gs.set[eglsRESETCOUNTERS] != 0)
2038 /* Reset all the counters related to performance over the run */
2039 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, runtime,
2040 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
2041 wcycle_set_reset_counters(wcycle, -1);
2042 if (!(cr->duty & DUTY_PME))
2044 /* Tell our PME node to reset its counters */
2045 gmx_pme_send_resetcounters(cr, step);
2047 /* Correct max_hours for the elapsed time */
2048 max_hours -= run_time/(60.0*60.0);
2049 bResetCountersHalfMaxH = FALSE;
2050 gs.set[eglsRESETCOUNTERS] = 0;
2054 /* End of main MD loop */
2058 runtime_end(runtime);
2060 if (bRerunMD && MASTER(cr))
2065 if (!(cr->duty & DUTY_PME))
2067 /* Tell the PME only node to finish */
2068 gmx_pme_send_finish(cr);
2073 if (ir->nstcalcenergy > 0 && !bRerunMD)
2075 print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
2076 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
2084 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
2086 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)));
2087 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
2090 if (pme_loadbal != NULL)
2092 pme_loadbal_done(pme_loadbal, cr, fplog,
2093 fr->nbv != NULL && fr->nbv->bUseGPU);
2096 if (shellfc && fplog)
2098 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
2099 (nconverged*100.0)/step_rel);
2100 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
2104 if (repl_ex_nst > 0 && MASTER(cr))
2106 print_replica_exchange_statistics(fplog, repl_ex);
2109 runtime->nsteps_done = step_rel;