2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
54 #include "md_support.h"
55 #include "md_logging.h"
68 #include "domdec_network.h"
69 #include "gromacs/gmxlib/topsort.h"
73 #include "compute_io.h"
74 #include "checkpoint.h"
75 #include "mtop_util.h"
76 #include "sighandler.h"
79 #include "pme_loadbal.h"
82 #include "types/nlistheuristics.h"
83 #include "types/iteratedconstraints.h"
84 #include "nbnxn_cuda_data_mgmt.h"
86 #include "gromacs/utility/gmxmpi.h"
87 #include "gromacs/fileio/confio.h"
88 #include "gromacs/fileio/trajectory_writing.h"
89 #include "gromacs/fileio/trnio.h"
90 #include "gromacs/fileio/trxio.h"
91 #include "gromacs/fileio/xtcio.h"
92 #include "gromacs/timing/wallcycle.h"
93 #include "gromacs/timing/walltime_accounting.h"
94 #include "gromacs/pulling/pull.h"
95 #include "gromacs/swap/swapcoords.h"
101 static void reset_all_counters(FILE *fplog, t_commrec *cr,
103 gmx_int64_t *step_rel, t_inputrec *ir,
104 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
105 gmx_walltime_accounting_t walltime_accounting,
106 nbnxn_cuda_ptr_t cu_nbv)
108 char sbuf[STEPSTRSIZE];
110 /* Reset all the counters related to performance over the run */
111 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
112 gmx_step_str(step, sbuf));
116 nbnxn_cuda_reset_timings(cu_nbv);
119 wallcycle_stop(wcycle, ewcRUN);
120 wallcycle_reset_all(wcycle);
121 if (DOMAINDECOMP(cr))
123 reset_dd_statistics_counters(cr->dd);
126 ir->init_step += *step_rel;
127 ir->nsteps -= *step_rel;
129 wallcycle_start(wcycle, ewcRUN);
130 walltime_accounting_start(walltime_accounting);
131 print_date_and_time(fplog, cr->nodeid, "Restarted time", walltime_accounting);
134 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
135 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
137 gmx_vsite_t *vsite, gmx_constr_t constr,
138 int stepout, t_inputrec *ir,
139 gmx_mtop_t *top_global,
141 t_state *state_global,
143 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
144 gmx_edsam_t ed, t_forcerec *fr,
145 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
146 real cpt_period, real max_hours,
147 const char gmx_unused *deviceOptions,
149 gmx_walltime_accounting_t walltime_accounting)
151 gmx_mdoutf_t outf = NULL;
152 gmx_int64_t step, step_rel;
154 double t, t0, lam0[efptNR];
155 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
156 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
157 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
158 bBornRadii, bStartingFromCpt;
159 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
160 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
161 bForceUpdate = FALSE, bCPT;
162 gmx_bool bMasterState;
163 int force_flags, cglo_flags;
164 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
169 t_state *bufstate = NULL;
170 matrix *scale_tot, pcoupl_mu, M, ebox;
173 gmx_repl_ex_t repl_ex = NULL;
176 t_mdebin *mdebin = NULL;
177 t_state *state = NULL;
178 rvec *f_global = NULL;
179 gmx_enerdata_t *enerd;
181 gmx_global_stat_t gstat;
182 gmx_update_t upd = NULL;
183 t_graph *graph = NULL;
185 gmx_groups_t *groups;
186 gmx_ekindata_t *ekind, *ekind_save;
187 gmx_shellfc_t shellfc;
188 int count, nconverged = 0;
191 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged, bNeedRepartition;
193 gmx_bool bResetCountersHalfMaxH = FALSE;
194 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
195 gmx_bool bUpdateDoLR;
199 real veta_save, scalevir, tracevir;
205 real saved_conserved_quantity = 0;
210 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
211 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
212 gmx_iterate_t iterate;
213 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
214 simulation stops. If equal to zero, don't
215 communicate any more between multisims.*/
216 /* PME load balancing data for GPU kernels */
217 pme_load_balancing_t pme_loadbal = NULL;
219 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
222 /* Temporary addition for FAHCORE checkpointing */
226 /* Check for special mdrun options */
227 bRerunMD = (Flags & MD_RERUN);
228 bAppend = (Flags & MD_APPENDFILES);
229 if (Flags & MD_RESETCOUNTERSHALFWAY)
233 /* Signal to reset the counters half the simulation steps. */
234 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
236 /* Signal to reset the counters halfway the simulation time. */
237 bResetCountersHalfMaxH = (max_hours > 0);
240 /* md-vv uses averaged full step velocities for T-control
241 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
242 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
244 if (bVV) /* to store the initial velocities while computing virial */
246 snew(cbuf, top_global->natoms);
248 /* all the iteratative cases - only if there are constraints */
249 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
250 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
251 false in this step. The correct value, true or false,
252 is set at each step, as it depends on the frequency of temperature
253 and pressure control.*/
254 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
258 /* Since we don't know if the frames read are related in any way,
259 * rebuild the neighborlist at every step.
262 ir->nstcalcenergy = 1;
266 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
268 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
269 bGStatEveryStep = (nstglobalcomm == 1);
271 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
274 "To reduce the energy communication with nstlist = -1\n"
275 "the neighbor list validity should not be checked at every step,\n"
276 "this means that exact integration is not guaranteed.\n"
277 "The neighbor list validity is checked after:\n"
278 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
279 "In most cases this will result in exact integration.\n"
280 "This reduces the energy communication by a factor of 2 to 3.\n"
281 "If you want less energy communication, set nstlist > 3.\n\n");
286 ir->nstxout_compressed = 0;
288 groups = &top_global->groups;
291 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
292 &(state_global->fep_state), lam0,
293 nrnb, top_global, &upd,
294 nfile, fnm, &outf, &mdebin,
295 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
297 clear_mat(total_vir);
299 /* Energy terms and groups */
301 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
303 if (DOMAINDECOMP(cr))
309 snew(f, top_global->natoms);
312 /* Kinetic energy data */
314 init_ekindata(fplog, top_global, &(ir->opts), ekind);
315 /* needed for iteration of constraints */
317 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
318 /* Copy the cos acceleration to the groups struct */
319 ekind->cosacc.cos_accel = ir->cos_accel;
321 gstat = global_stat_init(ir);
324 /* Check for polarizable models and flexible constraints */
325 shellfc = init_shell_flexcon(fplog, fr->cutoff_scheme == ecutsVERLET,
326 top_global, n_flexible_constraints(constr),
327 (ir->bContinuation ||
328 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
329 NULL : state_global->x);
333 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
334 set_deform_reference_box(upd,
335 deform_init_init_step_tpx,
336 deform_init_box_tpx);
337 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
341 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
342 if ((io > 2000) && MASTER(cr))
345 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
350 if (DOMAINDECOMP(cr))
352 top = dd_init_local_top(top_global);
355 dd_init_local_state(cr->dd, state_global, state);
357 if (DDMASTER(cr->dd) && ir->nstfout)
359 snew(f_global, state_global->natoms);
364 top = gmx_mtop_generate_local_top(top_global, ir);
366 forcerec_set_excl_load(fr, top);
368 state = serial_init_local_state(state_global);
371 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
375 set_vsite_top(vsite, top, mdatoms, cr);
378 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
380 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
385 make_local_shells(cr, mdatoms, shellfc);
388 setup_bonded_threading(fr, &top->idef);
391 if (DOMAINDECOMP(cr))
393 /* Distribute the charge groups over the nodes from the master node */
394 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
395 state_global, top_global, ir,
396 state, &f, mdatoms, top, fr,
397 vsite, shellfc, constr,
398 nrnb, wcycle, FALSE);
402 update_mdatoms(mdatoms, state->lambda[efptMASS]);
404 if (opt2bSet("-cpi", nfile, fnm))
406 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
410 bStateFromCP = FALSE;
415 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
422 /* Update mdebin with energy history if appending to output files */
423 if (Flags & MD_APPENDFILES)
425 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
429 /* We might have read an energy history from checkpoint,
430 * free the allocated memory and reset the counts.
432 done_energyhistory(&state_global->enerhist);
433 init_energyhistory(&state_global->enerhist);
436 /* Set the initial energy history in state by updating once */
437 update_energyhistory(&state_global->enerhist, mdebin);
440 /* Initialize constraints */
441 if (constr && !DOMAINDECOMP(cr))
443 set_constraints(constr, top, ir, mdatoms, cr);
448 /* We need to be sure replica exchange can only occur
449 * when the energies are current */
450 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
451 "repl_ex_nst", &repl_ex_nst);
452 /* This check needs to happen before inter-simulation
453 * signals are initialized, too */
455 if (repl_ex_nst > 0 && MASTER(cr))
457 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
458 repl_ex_nst, repl_ex_nex, repl_ex_seed);
461 /* PME tuning is only supported with GPUs or PME nodes and not with rerun or LJ-PME.
462 * With perturbed charges with soft-core we should not change the cut-off.
464 if ((Flags & MD_TUNEPME) &&
465 EEL_PME(fr->eeltype) &&
466 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
467 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
468 !bRerunMD && !EVDW_PME(fr->vdwtype))
470 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
472 if (cr->duty & DUTY_PME)
474 /* Start tuning right away, as we can't measure the load */
475 bPMETuneRunning = TRUE;
479 /* Separate PME nodes, we can measure the PP/PME load balance */
484 if (!ir->bContinuation && !bRerunMD)
486 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
488 /* Set the velocities of frozen particles to zero */
489 for (i = 0; i < mdatoms->homenr; i++)
491 for (m = 0; m < DIM; m++)
493 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
503 /* Constrain the initial coordinates and velocities */
504 do_constrain_first(fplog, constr, ir, mdatoms, state,
509 /* Construct the virtual sites for the initial configuration */
510 construct_vsites(vsite, state->x, ir->delta_t, NULL,
511 top->idef.iparams, top->idef.il,
512 fr->ePBC, fr->bMolPBC, cr, state->box);
518 /* set free energy calculation frequency as the minimum
519 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
520 nstfep = ir->fepvals->nstdhdl;
523 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
527 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
530 /* I'm assuming we need global communication the first time! MRS */
531 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
532 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
533 | (bVV ? CGLO_PRESSURE : 0)
534 | (bVV ? CGLO_CONSTRAINT : 0)
535 | (bRerunMD ? CGLO_RERUNMD : 0)
536 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
538 bSumEkinhOld = FALSE;
539 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
540 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
541 constr, NULL, FALSE, state->box,
542 top_global, &bSumEkinhOld, cglo_flags);
543 if (ir->eI == eiVVAK)
545 /* a second call to get the half step temperature initialized as well */
546 /* we do the same call as above, but turn the pressure off -- internally to
547 compute_globals, this is recognized as a velocity verlet half-step
548 kinetic energy calculation. This minimized excess variables, but
549 perhaps loses some logic?*/
551 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
552 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
553 constr, NULL, FALSE, state->box,
554 top_global, &bSumEkinhOld,
555 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
558 /* Calculate the initial half step temperature, and save the ekinh_old */
559 if (!(Flags & MD_STARTFROMCPT))
561 for (i = 0; (i < ir->opts.ngtc); i++)
563 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
568 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
569 and there is no previous step */
572 /* if using an iterative algorithm, we need to create a working directory for the state. */
575 bufstate = init_bufstate(state);
578 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
579 temperature control */
580 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
584 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
587 "RMS relative constraint deviation after constraining: %.2e\n",
588 constr_rmsd(constr, FALSE));
590 if (EI_STATE_VELOCITY(ir->eI))
592 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
596 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
597 " input trajectory '%s'\n\n",
598 *(top_global->name), opt2fn("-rerun", nfile, fnm));
601 fprintf(stderr, "Calculated time to finish depends on nsteps from "
602 "run input file,\nwhich may not correspond to the time "
603 "needed to process input trajectory.\n\n");
609 fprintf(stderr, "starting mdrun '%s'\n",
610 *(top_global->name));
613 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
617 sprintf(tbuf, "%s", "infinite");
619 if (ir->init_step > 0)
621 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
622 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
623 gmx_step_str(ir->init_step, sbuf2),
624 ir->init_step*ir->delta_t);
628 fprintf(stderr, "%s steps, %s ps.\n",
629 gmx_step_str(ir->nsteps, sbuf), tbuf);
632 fprintf(fplog, "\n");
635 walltime_accounting_start(walltime_accounting);
636 wallcycle_start(wcycle, ewcRUN);
637 print_start(fplog, cr, walltime_accounting, "mdrun");
639 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
641 chkpt_ret = fcCheckPointParallel( cr->nodeid,
645 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
650 /***********************************************************
654 ************************************************************/
656 /* if rerunMD then read coordinates and velocities from input trajectory */
659 if (getenv("GMX_FORCE_UPDATE"))
667 bNotLastFrame = read_first_frame(oenv, &status,
668 opt2fn("-rerun", nfile, fnm),
669 &rerun_fr, TRX_NEED_X | TRX_READ_V);
670 if (rerun_fr.natoms != top_global->natoms)
673 "Number of atoms in trajectory (%d) does not match the "
674 "run input file (%d)\n",
675 rerun_fr.natoms, top_global->natoms);
677 if (ir->ePBC != epbcNONE)
681 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);
683 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
685 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
692 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
695 if (ir->ePBC != epbcNONE)
697 /* Set the shift vectors.
698 * Necessary here when have a static box different from the tpr box.
700 calc_shifts(rerun_fr.box, fr->shift_vec);
704 /* loop over MD steps or if rerunMD to end of input trajectory */
706 /* Skip the first Nose-Hoover integration when we get the state from tpx */
707 bStateFromTPX = !bStateFromCP;
708 bInitStep = bFirstStep && (bStateFromTPX || bVV);
709 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
711 bSumEkinhOld = FALSE;
713 bNeedRepartition = FALSE;
715 init_global_signals(&gs, cr, ir, repl_ex_nst);
717 step = ir->init_step;
720 if (ir->nstlist == -1)
722 init_nlistheuristics(&nlh, bGStatEveryStep, step);
725 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
727 /* check how many steps are left in other sims */
728 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
732 /* and stop now if we should */
733 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
734 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
735 while (!bLastStep || (bRerunMD && bNotLastFrame))
738 wallcycle_start(wcycle, ewcSTEP);
744 step = rerun_fr.step;
745 step_rel = step - ir->init_step;
758 bLastStep = (step_rel == ir->nsteps);
759 t = t0 + step*ir->delta_t;
762 if (ir->efep != efepNO || ir->bSimTemp)
764 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
765 requiring different logic. */
767 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
768 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
769 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
770 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
771 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
776 update_annealing_target_temp(&(ir->opts), t);
781 if (!DOMAINDECOMP(cr) || MASTER(cr))
783 for (i = 0; i < state_global->natoms; i++)
785 copy_rvec(rerun_fr.x[i], state_global->x[i]);
789 for (i = 0; i < state_global->natoms; i++)
791 copy_rvec(rerun_fr.v[i], state_global->v[i]);
796 for (i = 0; i < state_global->natoms; i++)
798 clear_rvec(state_global->v[i]);
802 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
803 " Ekin, temperature and pressure are incorrect,\n"
804 " the virial will be incorrect when constraints are present.\n"
806 bRerunWarnNoV = FALSE;
810 copy_mat(rerun_fr.box, state_global->box);
811 copy_mat(state_global->box, state->box);
813 if (vsite && (Flags & MD_RERUN_VSITE))
815 if (DOMAINDECOMP(cr))
817 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
821 /* Following is necessary because the graph may get out of sync
822 * with the coordinates if we only have every N'th coordinate set
824 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
825 shift_self(graph, state->box, state->x);
827 construct_vsites(vsite, state->x, ir->delta_t, state->v,
828 top->idef.iparams, top->idef.il,
829 fr->ePBC, fr->bMolPBC, cr, state->box);
832 unshift_self(graph, state->box, state->x);
837 /* Stop Center of Mass motion */
838 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
842 /* for rerun MD always do Neighbour Searching */
843 bNS = (bFirstStep || ir->nstlist != 0);
848 /* Determine whether or not to do Neighbour Searching and LR */
849 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
851 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
852 (ir->nstlist == -1 && nlh.nabnsb > 0));
854 if (bNS && ir->nstlist == -1)
856 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
860 /* check whether we should stop because another simulation has
864 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
865 (multisim_nsteps != ir->nsteps) )
872 "Stopping simulation %d because another one has finished\n",
876 gs.sig[eglsCHKPT] = 1;
881 /* < 0 means stop at next step, > 0 means stop at next NS step */
882 if ( (gs.set[eglsSTOPCOND] < 0) ||
883 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
888 /* Determine whether or not to update the Born radii if doing GB */
889 bBornRadii = bFirstStep;
890 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
895 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
896 do_verbose = bVerbose &&
897 (step % stepout == 0 || bFirstStep || bLastStep);
899 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
907 bMasterState = FALSE;
908 /* Correct the new box if it is too skewed */
909 if (DYNAMIC_BOX(*ir))
911 if (correct_box(fplog, step, state->box, graph))
916 if (DOMAINDECOMP(cr) && bMasterState)
918 dd_collect_state(cr->dd, state, state_global);
922 if (DOMAINDECOMP(cr))
924 /* Repartition the domain decomposition */
925 wallcycle_start(wcycle, ewcDOMDEC);
926 dd_partition_system(fplog, step, cr,
927 bMasterState, nstglobalcomm,
928 state_global, top_global, ir,
929 state, &f, mdatoms, top, fr,
930 vsite, shellfc, constr,
932 do_verbose && !bPMETuneRunning);
933 wallcycle_stop(wcycle, ewcDOMDEC);
934 /* If using an iterative integrator, reallocate space to match the decomposition */
938 if (MASTER(cr) && do_log)
940 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
943 if (ir->efep != efepNO)
945 update_mdatoms(mdatoms, state->lambda[efptMASS]);
948 if ((bRerunMD && rerun_fr.bV) || bExchanged)
951 /* We need the kinetic energy at minus the half step for determining
952 * the full step kinetic energy and possibly for T-coupling.*/
953 /* This may not be quite working correctly yet . . . . */
954 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
955 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
956 constr, NULL, FALSE, state->box,
957 top_global, &bSumEkinhOld,
958 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
960 clear_mat(force_vir);
962 /* We write a checkpoint at this MD step when:
963 * either at an NS step when we signalled through gs,
964 * or at the last step (but not when we do not want confout),
965 * but never at the first step or with rerun.
967 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
968 (bLastStep && (Flags & MD_CONFOUT))) &&
969 step > ir->init_step && !bRerunMD);
972 gs.set[eglsCHKPT] = 0;
975 /* Determine the energy and pressure:
976 * at nstcalcenergy steps and at energy output steps (set below).
978 if (EI_VV(ir->eI) && (!bInitStep))
980 /* for vv, the first half of the integration actually corresponds
981 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
982 but the virial needs to be calculated on both the current step and the 'next' step. Future
983 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
985 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
986 bCalcVir = bCalcEner ||
987 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
991 bCalcEner = do_per_step(step, ir->nstcalcenergy);
992 bCalcVir = bCalcEner ||
993 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
996 /* Do we need global communication ? */
997 bGStat = (bCalcVir || bCalcEner || bStopCM ||
998 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
999 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1001 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1003 if (do_ene || do_log)
1010 /* these CGLO_ options remain the same throughout the iteration */
1011 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1012 (bGStat ? CGLO_GSTAT : 0)
1015 force_flags = (GMX_FORCE_STATECHANGED |
1016 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1017 GMX_FORCE_ALLFORCES |
1019 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1020 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1021 (bDoFEP ? GMX_FORCE_DHDL : 0)
1026 if (do_per_step(step, ir->nstcalclr))
1028 force_flags |= GMX_FORCE_DO_LR;
1034 /* Now is the time to relax the shells */
1035 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1036 ir, bNS, force_flags,
1039 state, f, force_vir, mdatoms,
1040 nrnb, wcycle, graph, groups,
1041 shellfc, fr, bBornRadii, t, mu_tot,
1043 mdoutf_get_fp_field(outf));
1053 /* The coordinates (x) are shifted (to get whole molecules)
1055 * This is parallellized as well, and does communication too.
1056 * Check comments in sim_util.c
1058 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1059 state->box, state->x, &state->hist,
1060 f, force_vir, mdatoms, enerd, fcd,
1061 state->lambda, graph,
1062 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1063 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1066 if (bVV && !bStartingFromCpt && !bRerunMD)
1067 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1069 if (ir->eI == eiVV && bInitStep)
1071 /* if using velocity verlet with full time step Ekin,
1072 * take the first half step only to compute the
1073 * virial for the first step. From there,
1074 * revert back to the initial coordinates
1075 * so that the input is actually the initial step.
1077 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1081 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1082 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1085 /* If we are using twin-range interactions where the long-range component
1086 * is only evaluated every nstcalclr>1 steps, we should do a special update
1087 * step to combine the long-range forces on these steps.
1088 * For nstcalclr=1 this is not done, since the forces would have been added
1089 * directly to the short-range forces already.
1091 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1093 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1094 f, bUpdateDoLR, fr->f_twin, fcd,
1095 ekind, M, upd, bInitStep, etrtVELOCITY1,
1096 cr, nrnb, constr, &top->idef);
1098 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1100 gmx_iterate_init(&iterate, TRUE);
1102 /* for iterations, we save these vectors, as we will be self-consistently iterating
1105 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1107 /* save the state */
1108 if (iterate.bIterationActive)
1110 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1113 bFirstIterate = TRUE;
1114 while (bFirstIterate || iterate.bIterationActive)
1116 if (iterate.bIterationActive)
1118 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1119 if (bFirstIterate && bTrotter)
1121 /* The first time through, we need a decent first estimate
1122 of veta(t+dt) to compute the constraints. Do
1123 this by computing the box volume part of the
1124 trotter integration at this time. Nothing else
1125 should be changed by this routine here. If
1126 !(first time), we start with the previous value
1129 veta_save = state->veta;
1130 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1131 vetanew = state->veta;
1132 state->veta = veta_save;
1137 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1139 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1140 state, fr->bMolPBC, graph, f,
1141 &top->idef, shake_vir,
1142 cr, nrnb, wcycle, upd, constr,
1143 TRUE, bCalcVir, vetanew);
1147 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1153 /* Need to unshift here if a do_force has been
1154 called in the previous step */
1155 unshift_self(graph, state->box, state->x);
1158 /* if VV, compute the pressure and constraints */
1159 /* For VV2, we strictly only need this if using pressure
1160 * control, but we really would like to have accurate pressures
1162 * Think about ways around this in the future?
1163 * For now, keep this choice in comments.
1165 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1166 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1168 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1169 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1171 bSumEkinhOld = TRUE;
1173 /* for vv, the first half of the integration actually corresponds to the previous step.
1174 So we need information from the last step in the first half of the integration */
1175 if (bGStat || do_per_step(step-1, nstglobalcomm))
1177 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1178 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1179 constr, NULL, FALSE, state->box,
1180 top_global, &bSumEkinhOld,
1183 | (bTemp ? CGLO_TEMPERATURE : 0)
1184 | (bPres ? CGLO_PRESSURE : 0)
1185 | (bPres ? CGLO_CONSTRAINT : 0)
1186 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1187 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1190 /* explanation of above:
1191 a) We compute Ekin at the full time step
1192 if 1) we are using the AveVel Ekin, and it's not the
1193 initial step, or 2) if we are using AveEkin, but need the full
1194 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1195 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1196 EkinAveVel because it's needed for the pressure */
1198 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1203 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1204 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1211 /* We need the kinetic energy at minus the half step for determining
1212 * the full step kinetic energy and possibly for T-coupling.*/
1213 /* This may not be quite working correctly yet . . . . */
1214 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1215 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1216 constr, NULL, FALSE, state->box,
1217 top_global, &bSumEkinhOld,
1218 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1223 if (iterate.bIterationActive &&
1224 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1225 state->veta, &vetanew))
1229 bFirstIterate = FALSE;
1232 if (bTrotter && !bInitStep)
1234 copy_mat(shake_vir, state->svir_prev);
1235 copy_mat(force_vir, state->fvir_prev);
1236 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1238 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1239 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1240 enerd->term[F_EKIN] = trace(ekind->ekin);
1243 /* if it's the initial step, we performed this first step just to get the constraint virial */
1244 if (bInitStep && ir->eI == eiVV)
1246 copy_rvecn(cbuf, state->v, 0, state->natoms);
1250 /* MRS -- now done iterating -- compute the conserved quantity */
1253 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1256 last_ekin = enerd->term[F_EKIN];
1258 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1260 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1262 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1265 sum_dhdl(enerd, state->lambda, ir->fepvals);
1269 /* ######## END FIRST UPDATE STEP ############## */
1270 /* ######## If doing VV, we now have v(dt) ###### */
1273 /* perform extended ensemble sampling in lambda - we don't
1274 actually move to the new state before outputting
1275 statistics, but if performing simulated tempering, we
1276 do update the velocities and the tau_t. */
1278 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1279 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1280 copy_df_history(&state_global->dfhist, &state->dfhist);
1283 /* Now we have the energies and forces corresponding to the
1284 * coordinates at time t. We must output all of this before
1287 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1288 ir, state, state_global, top_global, fr,
1289 outf, mdebin, ekind, f, f_global,
1291 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1294 /* kludge -- virial is lost with restart for NPT control. Must restart */
1295 if (bStartingFromCpt && bVV)
1297 copy_mat(state->svir_prev, shake_vir);
1298 copy_mat(state->fvir_prev, force_vir);
1301 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1303 /* Check whether everything is still allright */
1304 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1305 #ifdef GMX_THREAD_MPI
1310 /* this is just make gs.sig compatible with the hack
1311 of sending signals around by MPI_Reduce with together with
1313 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1315 gs.sig[eglsSTOPCOND] = 1;
1317 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1319 gs.sig[eglsSTOPCOND] = -1;
1321 /* < 0 means stop at next step, > 0 means stop at next NS step */
1325 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1326 gmx_get_signal_name(),
1327 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1331 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1332 gmx_get_signal_name(),
1333 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1335 handled_stop_condition = (int)gmx_get_stop_condition();
1337 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1338 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1339 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1341 /* Signal to terminate the run */
1342 gs.sig[eglsSTOPCOND] = 1;
1345 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1347 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1350 if (bResetCountersHalfMaxH && MASTER(cr) &&
1351 elapsed_time > max_hours*60.0*60.0*0.495)
1353 gs.sig[eglsRESETCOUNTERS] = 1;
1356 if (ir->nstlist == -1 && !bRerunMD)
1358 /* When bGStatEveryStep=FALSE, global_stat is only called
1359 * when we check the atom displacements, not at NS steps.
1360 * This means that also the bonded interaction count check is not
1361 * performed immediately after NS. Therefore a few MD steps could
1362 * be performed with missing interactions.
1363 * But wrong energies are never written to file,
1364 * since energies are only written after global_stat
1367 if (step >= nlh.step_nscheck)
1369 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1370 nlh.scale_tot, state->x);
1374 /* This is not necessarily true,
1375 * but step_nscheck is determined quite conservatively.
1381 /* In parallel we only have to check for checkpointing in steps
1382 * where we do global communication,
1383 * otherwise the other nodes don't know.
1385 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1388 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1389 gs.set[eglsCHKPT] == 0)
1391 gs.sig[eglsCHKPT] = 1;
1394 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1399 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1401 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1403 gmx_bool bIfRandomize;
1404 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1405 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1406 if (constr && bIfRandomize)
1408 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1409 state, fr->bMolPBC, graph, f,
1410 &top->idef, tmp_vir,
1411 cr, nrnb, wcycle, upd, constr,
1412 TRUE, bCalcVir, vetanew);
1417 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1419 gmx_iterate_init(&iterate, TRUE);
1420 /* for iterations, we save these vectors, as we will be redoing the calculations */
1421 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1424 bFirstIterate = TRUE;
1425 while (bFirstIterate || iterate.bIterationActive)
1427 /* We now restore these vectors to redo the calculation with improved extended variables */
1428 if (iterate.bIterationActive)
1430 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1433 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1434 so scroll down for that logic */
1436 /* ######### START SECOND UPDATE STEP ################# */
1437 /* Box is changed in update() when we do pressure coupling,
1438 * but we should still use the old box for energy corrections and when
1439 * writing it to the energy file, so it matches the trajectory files for
1440 * the same timestep above. Make a copy in a separate array.
1442 copy_mat(state->box, lastbox);
1447 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1449 wallcycle_start(wcycle, ewcUPDATE);
1450 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1453 if (iterate.bIterationActive)
1461 /* we use a new value of scalevir to converge the iterations faster */
1462 scalevir = tracevir/trace(shake_vir);
1464 msmul(shake_vir, scalevir, shake_vir);
1465 m_add(force_vir, shake_vir, total_vir);
1466 clear_mat(shake_vir);
1468 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1469 /* We can only do Berendsen coupling after we have summed
1470 * the kinetic energy or virial. Since the happens
1471 * in global_state after update, we should only do it at
1472 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1477 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1478 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1483 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1485 /* velocity half-step update */
1486 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1487 bUpdateDoLR, fr->f_twin, fcd,
1488 ekind, M, upd, FALSE, etrtVELOCITY2,
1489 cr, nrnb, constr, &top->idef);
1492 /* Above, initialize just copies ekinh into ekin,
1493 * it doesn't copy position (for VV),
1494 * and entire integrator for MD.
1497 if (ir->eI == eiVVAK)
1499 copy_rvecn(state->x, cbuf, 0, state->natoms);
1501 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1503 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1504 bUpdateDoLR, fr->f_twin, fcd,
1505 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1506 wallcycle_stop(wcycle, ewcUPDATE);
1508 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1509 fr->bMolPBC, graph, f,
1510 &top->idef, shake_vir,
1511 cr, nrnb, wcycle, upd, constr,
1512 FALSE, bCalcVir, state->veta);
1514 if (ir->eI == eiVVAK)
1516 /* erase F_EKIN and F_TEMP here? */
1517 /* just compute the kinetic energy at the half step to perform a trotter step */
1518 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1519 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1520 constr, NULL, FALSE, lastbox,
1521 top_global, &bSumEkinhOld,
1522 cglo_flags | CGLO_TEMPERATURE
1524 wallcycle_start(wcycle, ewcUPDATE);
1525 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1526 /* now we know the scaling, we can compute the positions again again */
1527 copy_rvecn(cbuf, state->x, 0, state->natoms);
1529 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1531 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1532 bUpdateDoLR, fr->f_twin, fcd,
1533 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1534 wallcycle_stop(wcycle, ewcUPDATE);
1536 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1537 /* are the small terms in the shake_vir here due
1538 * to numerical errors, or are they important
1539 * physically? I'm thinking they are just errors, but not completely sure.
1540 * For now, will call without actually constraining, constr=NULL*/
1541 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1542 state, fr->bMolPBC, graph, f,
1543 &top->idef, tmp_vir,
1544 cr, nrnb, wcycle, upd, NULL,
1550 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1553 if (fr->bSepDVDL && fplog && do_log)
1555 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1559 /* this factor or 2 correction is necessary
1560 because half of the constraint force is removed
1561 in the vv step, so we have to double it. See
1562 the Redmine issue #1255. It is not yet clear
1563 if the factor of 2 is exact, or just a very
1564 good approximation, and this will be
1565 investigated. The next step is to see if this
1566 can be done adding a dhdl contribution from the
1567 rattle step, but this is somewhat more
1568 complicated with the current code. Will be
1569 investigated, hopefully for 4.6.3. However,
1570 this current solution is much better than
1571 having it completely wrong.
1573 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1577 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1582 /* Need to unshift here */
1583 unshift_self(graph, state->box, state->x);
1588 wallcycle_start(wcycle, ewcVSITECONSTR);
1591 shift_self(graph, state->box, state->x);
1593 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1594 top->idef.iparams, top->idef.il,
1595 fr->ePBC, fr->bMolPBC, cr, state->box);
1599 unshift_self(graph, state->box, state->x);
1601 wallcycle_stop(wcycle, ewcVSITECONSTR);
1604 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1605 /* With Leap-Frog we can skip compute_globals at
1606 * non-communication steps, but we need to calculate
1607 * the kinetic energy one step before communication.
1609 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1611 if (ir->nstlist == -1 && bFirstIterate)
1613 gs.sig[eglsNABNSB] = nlh.nabnsb;
1615 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1616 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1618 bFirstIterate ? &gs : NULL,
1619 (step_rel % gs.nstms == 0) &&
1620 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1622 top_global, &bSumEkinhOld,
1624 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1625 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1626 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1627 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1628 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1629 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1632 if (ir->nstlist == -1 && bFirstIterate)
1634 nlh.nabnsb = gs.set[eglsNABNSB];
1635 gs.set[eglsNABNSB] = 0;
1638 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1639 /* ############# END CALC EKIN AND PRESSURE ################# */
1641 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1642 the virial that should probably be addressed eventually. state->veta has better properies,
1643 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1644 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1646 if (iterate.bIterationActive &&
1647 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1648 trace(shake_vir), &tracevir))
1652 bFirstIterate = FALSE;
1655 if (!bVV || bRerunMD)
1657 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1658 sum_dhdl(enerd, state->lambda, ir->fepvals);
1660 update_box(fplog, step, ir, mdatoms, state, f,
1661 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1663 /* ################# END UPDATE STEP 2 ################# */
1664 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1666 /* The coordinates (x) were unshifted in update */
1669 /* We will not sum ekinh_old,
1670 * so signal that we still have to do it.
1672 bSumEkinhOld = TRUE;
1675 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1677 /* use the directly determined last velocity, not actually the averaged half steps */
1678 if (bTrotter && ir->eI == eiVV)
1680 enerd->term[F_EKIN] = last_ekin;
1682 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1686 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1690 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1692 /* ######### END PREPARING EDR OUTPUT ########### */
1697 gmx_bool do_dr, do_or;
1699 if (fplog && do_log && bDoExpanded)
1701 /* only needed if doing expanded ensemble */
1702 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1703 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1705 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1709 upd_mdebin(mdebin, bDoDHDL, TRUE,
1710 t, mdatoms->tmass, enerd, state,
1711 ir->fepvals, ir->expandedvals, lastbox,
1712 shake_vir, force_vir, total_vir, pres,
1713 ekind, mu_tot, constr);
1717 upd_mdebin_step(mdebin);
1720 do_dr = do_per_step(step, ir->nstdisreout);
1721 do_or = do_per_step(step, ir->nstorireout);
1723 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1725 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1727 if (ir->ePull != epullNO)
1729 pull_print_output(ir->pull, step, t);
1732 if (do_per_step(step, ir->nstlog))
1734 if (fflush(fplog) != 0)
1736 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1742 /* Have to do this part _after_ outputting the logfile and the edr file */
1743 /* Gets written into the state at the beginning of next loop*/
1744 state->fep_state = lamnew;
1746 /* Print the remaining wall clock time for the run */
1747 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1751 fprintf(stderr, "\n");
1753 print_time(stderr, walltime_accounting, step, ir, cr);
1756 /* Ion/water position swapping.
1757 * Not done in last step since trajectory writing happens before this call
1758 * in the MD loop and exchanges would be lost anyway. */
1759 bNeedRepartition = FALSE;
1760 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1761 do_per_step(step, ir->swap->nstswap))
1763 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1764 bRerunMD ? rerun_fr.x : state->x,
1765 bRerunMD ? rerun_fr.box : state->box,
1766 top_global, MASTER(cr) && bVerbose, bRerunMD);
1768 if (bNeedRepartition && DOMAINDECOMP(cr))
1770 dd_collect_state(cr->dd, state, state_global);
1774 /* Replica exchange */
1776 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1777 do_per_step(step, repl_ex_nst))
1779 bExchanged = replica_exchange(fplog, cr, repl_ex,
1780 state_global, enerd,
1784 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1786 dd_partition_system(fplog, step, cr, TRUE, 1,
1787 state_global, top_global, ir,
1788 state, &f, mdatoms, top, fr,
1789 vsite, shellfc, constr,
1790 nrnb, wcycle, FALSE);
1795 bStartingFromCpt = FALSE;
1797 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1798 /* With all integrators, except VV, we need to retain the pressure
1799 * at the current step for coupling at the next step.
1801 if ((state->flags & (1<<estPRES_PREV)) &&
1803 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1805 /* Store the pressure in t_state for pressure coupling
1806 * at the next MD step.
1808 copy_mat(pres, state->pres_prev);
1811 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1813 if ( (membed != NULL) && (!bLastStep) )
1815 rescale_membed(step_rel, membed, state_global->x);
1822 /* read next frame from input trajectory */
1823 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1828 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1832 if (!bRerunMD || !rerun_fr.bStep)
1834 /* increase the MD step number */
1839 cycles = wallcycle_stop(wcycle, ewcSTEP);
1840 if (DOMAINDECOMP(cr) && wcycle)
1842 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1845 if (bPMETuneRunning || bPMETuneTry)
1847 /* PME grid + cut-off optimization with GPUs or PME nodes */
1849 /* Count the total cycles over the last steps */
1850 cycles_pmes += cycles;
1852 /* We can only switch cut-off at NS steps */
1853 if (step % ir->nstlist == 0)
1855 /* PME grid + cut-off optimization with GPUs or PME nodes */
1858 if (DDMASTER(cr->dd))
1860 /* PME node load is too high, start tuning */
1861 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1863 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1865 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1867 bPMETuneTry = FALSE;
1870 if (bPMETuneRunning)
1872 /* init_step might not be a multiple of nstlist,
1873 * but the first cycle is always skipped anyhow.
1876 pme_load_balance(pme_loadbal, cr,
1877 (bVerbose && MASTER(cr)) ? stderr : NULL,
1879 ir, state, cycles_pmes,
1880 fr->ic, fr->nbv, &fr->pmedata,
1883 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1884 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1885 fr->rlist = fr->ic->rlist;
1886 fr->rlistlong = fr->ic->rlistlong;
1887 fr->rcoulomb = fr->ic->rcoulomb;
1888 fr->rvdw = fr->ic->rvdw;
1894 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1895 gs.set[eglsRESETCOUNTERS] != 0)
1897 /* Reset all the counters related to performance over the run */
1898 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1899 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1900 wcycle_set_reset_counters(wcycle, -1);
1901 if (!(cr->duty & DUTY_PME))
1903 /* Tell our PME node to reset its counters */
1904 gmx_pme_send_resetcounters(cr, step);
1906 /* Correct max_hours for the elapsed time */
1907 max_hours -= elapsed_time/(60.0*60.0);
1908 bResetCountersHalfMaxH = FALSE;
1909 gs.set[eglsRESETCOUNTERS] = 0;
1913 /* End of main MD loop */
1916 /* Stop measuring walltime */
1917 walltime_accounting_end(walltime_accounting);
1919 if (bRerunMD && MASTER(cr))
1924 if (!(cr->duty & DUTY_PME))
1926 /* Tell the PME only node to finish */
1927 gmx_pme_send_finish(cr);
1932 if (ir->nstcalcenergy > 0 && !bRerunMD)
1934 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1935 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1942 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1944 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)));
1945 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1948 if (pme_loadbal != NULL)
1950 pme_loadbal_done(pme_loadbal, cr, fplog,
1951 fr->nbv != NULL && fr->nbv->bUseGPU);
1954 if (shellfc && fplog)
1956 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1957 (nconverged*100.0)/step_rel);
1958 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1962 if (repl_ex_nst > 0 && MASTER(cr))
1964 print_replica_exchange_statistics(fplog, repl_ex);
1967 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);