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.
462 * PME tuning is not supported with PME only for LJ and not for Coulomb.
464 if ((Flags & MD_TUNEPME) &&
465 EEL_PME(fr->eeltype) &&
466 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
469 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
471 if (cr->duty & DUTY_PME)
473 /* Start tuning right away, as we can't measure the load */
474 bPMETuneRunning = TRUE;
478 /* Separate PME nodes, we can measure the PP/PME load balance */
483 if (!ir->bContinuation && !bRerunMD)
485 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
487 /* Set the velocities of frozen particles to zero */
488 for (i = 0; i < mdatoms->homenr; i++)
490 for (m = 0; m < DIM; m++)
492 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
502 /* Constrain the initial coordinates and velocities */
503 do_constrain_first(fplog, constr, ir, mdatoms, state,
508 /* Construct the virtual sites for the initial configuration */
509 construct_vsites(vsite, state->x, ir->delta_t, NULL,
510 top->idef.iparams, top->idef.il,
511 fr->ePBC, fr->bMolPBC, cr, state->box);
517 /* set free energy calculation frequency as the minimum
518 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
519 nstfep = ir->fepvals->nstdhdl;
522 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
526 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
529 /* I'm assuming we need global communication the first time! MRS */
530 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
531 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
532 | (bVV ? CGLO_PRESSURE : 0)
533 | (bVV ? CGLO_CONSTRAINT : 0)
534 | (bRerunMD ? CGLO_RERUNMD : 0)
535 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
537 bSumEkinhOld = FALSE;
538 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
539 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
540 constr, NULL, FALSE, state->box,
541 top_global, &bSumEkinhOld, cglo_flags);
542 if (ir->eI == eiVVAK)
544 /* a second call to get the half step temperature initialized as well */
545 /* we do the same call as above, but turn the pressure off -- internally to
546 compute_globals, this is recognized as a velocity verlet half-step
547 kinetic energy calculation. This minimized excess variables, but
548 perhaps loses some logic?*/
550 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
551 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
552 constr, NULL, FALSE, state->box,
553 top_global, &bSumEkinhOld,
554 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
557 /* Calculate the initial half step temperature, and save the ekinh_old */
558 if (!(Flags & MD_STARTFROMCPT))
560 for (i = 0; (i < ir->opts.ngtc); i++)
562 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
567 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
568 and there is no previous step */
571 /* if using an iterative algorithm, we need to create a working directory for the state. */
574 bufstate = init_bufstate(state);
577 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
578 temperature control */
579 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
583 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
586 "RMS relative constraint deviation after constraining: %.2e\n",
587 constr_rmsd(constr, FALSE));
589 if (EI_STATE_VELOCITY(ir->eI))
591 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
595 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
596 " input trajectory '%s'\n\n",
597 *(top_global->name), opt2fn("-rerun", nfile, fnm));
600 fprintf(stderr, "Calculated time to finish depends on nsteps from "
601 "run input file,\nwhich may not correspond to the time "
602 "needed to process input trajectory.\n\n");
608 fprintf(stderr, "starting mdrun '%s'\n",
609 *(top_global->name));
612 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
616 sprintf(tbuf, "%s", "infinite");
618 if (ir->init_step > 0)
620 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
621 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
622 gmx_step_str(ir->init_step, sbuf2),
623 ir->init_step*ir->delta_t);
627 fprintf(stderr, "%s steps, %s ps.\n",
628 gmx_step_str(ir->nsteps, sbuf), tbuf);
631 fprintf(fplog, "\n");
634 walltime_accounting_start(walltime_accounting);
635 wallcycle_start(wcycle, ewcRUN);
636 print_start(fplog, cr, walltime_accounting, "mdrun");
638 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
640 chkpt_ret = fcCheckPointParallel( cr->nodeid,
644 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
649 /***********************************************************
653 ************************************************************/
655 /* if rerunMD then read coordinates and velocities from input trajectory */
658 if (getenv("GMX_FORCE_UPDATE"))
666 bNotLastFrame = read_first_frame(oenv, &status,
667 opt2fn("-rerun", nfile, fnm),
668 &rerun_fr, TRX_NEED_X | TRX_READ_V);
669 if (rerun_fr.natoms != top_global->natoms)
672 "Number of atoms in trajectory (%d) does not match the "
673 "run input file (%d)\n",
674 rerun_fr.natoms, top_global->natoms);
676 if (ir->ePBC != epbcNONE)
680 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);
682 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
684 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
691 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
694 if (ir->ePBC != epbcNONE)
696 /* Set the shift vectors.
697 * Necessary here when have a static box different from the tpr box.
699 calc_shifts(rerun_fr.box, fr->shift_vec);
703 /* loop over MD steps or if rerunMD to end of input trajectory */
705 /* Skip the first Nose-Hoover integration when we get the state from tpx */
706 bStateFromTPX = !bStateFromCP;
707 bInitStep = bFirstStep && (bStateFromTPX || bVV);
708 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
710 bSumEkinhOld = FALSE;
712 bNeedRepartition = FALSE;
714 init_global_signals(&gs, cr, ir, repl_ex_nst);
716 step = ir->init_step;
719 if (ir->nstlist == -1)
721 init_nlistheuristics(&nlh, bGStatEveryStep, step);
724 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
726 /* check how many steps are left in other sims */
727 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
731 /* and stop now if we should */
732 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
733 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
734 while (!bLastStep || (bRerunMD && bNotLastFrame))
737 wallcycle_start(wcycle, ewcSTEP);
743 step = rerun_fr.step;
744 step_rel = step - ir->init_step;
757 bLastStep = (step_rel == ir->nsteps);
758 t = t0 + step*ir->delta_t;
761 if (ir->efep != efepNO || ir->bSimTemp)
763 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
764 requiring different logic. */
766 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
767 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
768 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
769 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
770 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
775 update_annealing_target_temp(&(ir->opts), t);
780 if (!DOMAINDECOMP(cr) || MASTER(cr))
782 for (i = 0; i < state_global->natoms; i++)
784 copy_rvec(rerun_fr.x[i], state_global->x[i]);
788 for (i = 0; i < state_global->natoms; i++)
790 copy_rvec(rerun_fr.v[i], state_global->v[i]);
795 for (i = 0; i < state_global->natoms; i++)
797 clear_rvec(state_global->v[i]);
801 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
802 " Ekin, temperature and pressure are incorrect,\n"
803 " the virial will be incorrect when constraints are present.\n"
805 bRerunWarnNoV = FALSE;
809 copy_mat(rerun_fr.box, state_global->box);
810 copy_mat(state_global->box, state->box);
812 if (vsite && (Flags & MD_RERUN_VSITE))
814 if (DOMAINDECOMP(cr))
816 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
820 /* Following is necessary because the graph may get out of sync
821 * with the coordinates if we only have every N'th coordinate set
823 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
824 shift_self(graph, state->box, state->x);
826 construct_vsites(vsite, state->x, ir->delta_t, state->v,
827 top->idef.iparams, top->idef.il,
828 fr->ePBC, fr->bMolPBC, cr, state->box);
831 unshift_self(graph, state->box, state->x);
836 /* Stop Center of Mass motion */
837 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
841 /* for rerun MD always do Neighbour Searching */
842 bNS = (bFirstStep || ir->nstlist != 0);
847 /* Determine whether or not to do Neighbour Searching and LR */
848 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
850 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
851 (ir->nstlist == -1 && nlh.nabnsb > 0));
853 if (bNS && ir->nstlist == -1)
855 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
859 /* check whether we should stop because another simulation has
863 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
864 (multisim_nsteps != ir->nsteps) )
871 "Stopping simulation %d because another one has finished\n",
875 gs.sig[eglsCHKPT] = 1;
880 /* < 0 means stop at next step, > 0 means stop at next NS step */
881 if ( (gs.set[eglsSTOPCOND] < 0) ||
882 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
887 /* Determine whether or not to update the Born radii if doing GB */
888 bBornRadii = bFirstStep;
889 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
894 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
895 do_verbose = bVerbose &&
896 (step % stepout == 0 || bFirstStep || bLastStep);
898 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
906 bMasterState = FALSE;
907 /* Correct the new box if it is too skewed */
908 if (DYNAMIC_BOX(*ir))
910 if (correct_box(fplog, step, state->box, graph))
915 if (DOMAINDECOMP(cr) && bMasterState)
917 dd_collect_state(cr->dd, state, state_global);
921 if (DOMAINDECOMP(cr))
923 /* Repartition the domain decomposition */
924 wallcycle_start(wcycle, ewcDOMDEC);
925 dd_partition_system(fplog, step, cr,
926 bMasterState, nstglobalcomm,
927 state_global, top_global, ir,
928 state, &f, mdatoms, top, fr,
929 vsite, shellfc, constr,
931 do_verbose && !bPMETuneRunning);
932 wallcycle_stop(wcycle, ewcDOMDEC);
933 /* If using an iterative integrator, reallocate space to match the decomposition */
937 if (MASTER(cr) && do_log)
939 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
942 if (ir->efep != efepNO)
944 update_mdatoms(mdatoms, state->lambda[efptMASS]);
947 if ((bRerunMD && rerun_fr.bV) || bExchanged)
950 /* We need the kinetic energy at minus the half step for determining
951 * the full step kinetic energy and possibly for T-coupling.*/
952 /* This may not be quite working correctly yet . . . . */
953 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
954 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
955 constr, NULL, FALSE, state->box,
956 top_global, &bSumEkinhOld,
957 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
959 clear_mat(force_vir);
961 /* We write a checkpoint at this MD step when:
962 * either at an NS step when we signalled through gs,
963 * or at the last step (but not when we do not want confout),
964 * but never at the first step or with rerun.
966 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
967 (bLastStep && (Flags & MD_CONFOUT))) &&
968 step > ir->init_step && !bRerunMD);
971 gs.set[eglsCHKPT] = 0;
974 /* Determine the energy and pressure:
975 * at nstcalcenergy steps and at energy output steps (set below).
977 if (EI_VV(ir->eI) && (!bInitStep))
979 /* for vv, the first half of the integration actually corresponds
980 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
981 but the virial needs to be calculated on both the current step and the 'next' step. Future
982 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
984 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
985 bCalcVir = bCalcEner ||
986 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
990 bCalcEner = do_per_step(step, ir->nstcalcenergy);
991 bCalcVir = bCalcEner ||
992 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
995 /* Do we need global communication ? */
996 bGStat = (bCalcVir || bCalcEner || bStopCM ||
997 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
998 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1000 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1002 if (do_ene || do_log)
1009 /* these CGLO_ options remain the same throughout the iteration */
1010 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1011 (bGStat ? CGLO_GSTAT : 0)
1014 force_flags = (GMX_FORCE_STATECHANGED |
1015 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1016 GMX_FORCE_ALLFORCES |
1018 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1019 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1020 (bDoFEP ? GMX_FORCE_DHDL : 0)
1025 if (do_per_step(step, ir->nstcalclr))
1027 force_flags |= GMX_FORCE_DO_LR;
1033 /* Now is the time to relax the shells */
1034 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1035 ir, bNS, force_flags,
1038 state, f, force_vir, mdatoms,
1039 nrnb, wcycle, graph, groups,
1040 shellfc, fr, bBornRadii, t, mu_tot,
1042 mdoutf_get_fp_field(outf));
1052 /* The coordinates (x) are shifted (to get whole molecules)
1054 * This is parallellized as well, and does communication too.
1055 * Check comments in sim_util.c
1057 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1058 state->box, state->x, &state->hist,
1059 f, force_vir, mdatoms, enerd, fcd,
1060 state->lambda, graph,
1061 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1062 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1065 if (bVV && !bStartingFromCpt && !bRerunMD)
1066 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1068 if (ir->eI == eiVV && bInitStep)
1070 /* if using velocity verlet with full time step Ekin,
1071 * take the first half step only to compute the
1072 * virial for the first step. From there,
1073 * revert back to the initial coordinates
1074 * so that the input is actually the initial step.
1076 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1080 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1081 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1084 /* If we are using twin-range interactions where the long-range component
1085 * is only evaluated every nstcalclr>1 steps, we should do a special update
1086 * step to combine the long-range forces on these steps.
1087 * For nstcalclr=1 this is not done, since the forces would have been added
1088 * directly to the short-range forces already.
1090 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1092 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1093 f, bUpdateDoLR, fr->f_twin, fcd,
1094 ekind, M, upd, bInitStep, etrtVELOCITY1,
1095 cr, nrnb, constr, &top->idef);
1097 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1099 gmx_iterate_init(&iterate, TRUE);
1101 /* for iterations, we save these vectors, as we will be self-consistently iterating
1104 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1106 /* save the state */
1107 if (iterate.bIterationActive)
1109 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1112 bFirstIterate = TRUE;
1113 while (bFirstIterate || iterate.bIterationActive)
1115 if (iterate.bIterationActive)
1117 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1118 if (bFirstIterate && bTrotter)
1120 /* The first time through, we need a decent first estimate
1121 of veta(t+dt) to compute the constraints. Do
1122 this by computing the box volume part of the
1123 trotter integration at this time. Nothing else
1124 should be changed by this routine here. If
1125 !(first time), we start with the previous value
1128 veta_save = state->veta;
1129 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1130 vetanew = state->veta;
1131 state->veta = veta_save;
1136 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1138 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1139 state, fr->bMolPBC, graph, f,
1140 &top->idef, shake_vir,
1141 cr, nrnb, wcycle, upd, constr,
1142 TRUE, bCalcVir, vetanew);
1146 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1152 /* Need to unshift here if a do_force has been
1153 called in the previous step */
1154 unshift_self(graph, state->box, state->x);
1157 /* if VV, compute the pressure and constraints */
1158 /* For VV2, we strictly only need this if using pressure
1159 * control, but we really would like to have accurate pressures
1161 * Think about ways around this in the future?
1162 * For now, keep this choice in comments.
1164 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1165 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1167 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1168 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1170 bSumEkinhOld = TRUE;
1172 /* for vv, the first half of the integration actually corresponds to the previous step.
1173 So we need information from the last step in the first half of the integration */
1174 if (bGStat || do_per_step(step-1, nstglobalcomm))
1176 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1177 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1178 constr, NULL, FALSE, state->box,
1179 top_global, &bSumEkinhOld,
1182 | (bTemp ? CGLO_TEMPERATURE : 0)
1183 | (bPres ? CGLO_PRESSURE : 0)
1184 | (bPres ? CGLO_CONSTRAINT : 0)
1185 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1186 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1189 /* explanation of above:
1190 a) We compute Ekin at the full time step
1191 if 1) we are using the AveVel Ekin, and it's not the
1192 initial step, or 2) if we are using AveEkin, but need the full
1193 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1194 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1195 EkinAveVel because it's needed for the pressure */
1197 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1202 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1203 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1210 /* We need the kinetic energy at minus the half step for determining
1211 * the full step kinetic energy and possibly for T-coupling.*/
1212 /* This may not be quite working correctly yet . . . . */
1213 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1214 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1215 constr, NULL, FALSE, state->box,
1216 top_global, &bSumEkinhOld,
1217 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1222 if (iterate.bIterationActive &&
1223 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1224 state->veta, &vetanew))
1228 bFirstIterate = FALSE;
1231 if (bTrotter && !bInitStep)
1233 copy_mat(shake_vir, state->svir_prev);
1234 copy_mat(force_vir, state->fvir_prev);
1235 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1237 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1238 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1239 enerd->term[F_EKIN] = trace(ekind->ekin);
1242 /* if it's the initial step, we performed this first step just to get the constraint virial */
1243 if (bInitStep && ir->eI == eiVV)
1245 copy_rvecn(cbuf, state->v, 0, state->natoms);
1249 /* MRS -- now done iterating -- compute the conserved quantity */
1252 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1255 last_ekin = enerd->term[F_EKIN];
1257 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1259 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1261 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1264 sum_dhdl(enerd, state->lambda, ir->fepvals);
1268 /* ######## END FIRST UPDATE STEP ############## */
1269 /* ######## If doing VV, we now have v(dt) ###### */
1272 /* perform extended ensemble sampling in lambda - we don't
1273 actually move to the new state before outputting
1274 statistics, but if performing simulated tempering, we
1275 do update the velocities and the tau_t. */
1277 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1278 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1279 copy_df_history(&state_global->dfhist, &state->dfhist);
1282 /* Now we have the energies and forces corresponding to the
1283 * coordinates at time t. We must output all of this before
1286 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1287 ir, state, state_global, top_global, fr,
1288 outf, mdebin, ekind, f, f_global,
1290 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1293 /* kludge -- virial is lost with restart for NPT control. Must restart */
1294 if (bStartingFromCpt && bVV)
1296 copy_mat(state->svir_prev, shake_vir);
1297 copy_mat(state->fvir_prev, force_vir);
1300 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1302 /* Check whether everything is still allright */
1303 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1304 #ifdef GMX_THREAD_MPI
1309 /* this is just make gs.sig compatible with the hack
1310 of sending signals around by MPI_Reduce with together with
1312 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1314 gs.sig[eglsSTOPCOND] = 1;
1316 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1318 gs.sig[eglsSTOPCOND] = -1;
1320 /* < 0 means stop at next step, > 0 means stop at next NS step */
1324 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1325 gmx_get_signal_name(),
1326 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1330 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1331 gmx_get_signal_name(),
1332 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1334 handled_stop_condition = (int)gmx_get_stop_condition();
1336 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1337 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1338 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1340 /* Signal to terminate the run */
1341 gs.sig[eglsSTOPCOND] = 1;
1344 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1346 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1349 if (bResetCountersHalfMaxH && MASTER(cr) &&
1350 elapsed_time > max_hours*60.0*60.0*0.495)
1352 gs.sig[eglsRESETCOUNTERS] = 1;
1355 if (ir->nstlist == -1 && !bRerunMD)
1357 /* When bGStatEveryStep=FALSE, global_stat is only called
1358 * when we check the atom displacements, not at NS steps.
1359 * This means that also the bonded interaction count check is not
1360 * performed immediately after NS. Therefore a few MD steps could
1361 * be performed with missing interactions.
1362 * But wrong energies are never written to file,
1363 * since energies are only written after global_stat
1366 if (step >= nlh.step_nscheck)
1368 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1369 nlh.scale_tot, state->x);
1373 /* This is not necessarily true,
1374 * but step_nscheck is determined quite conservatively.
1380 /* In parallel we only have to check for checkpointing in steps
1381 * where we do global communication,
1382 * otherwise the other nodes don't know.
1384 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1387 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1388 gs.set[eglsCHKPT] == 0)
1390 gs.sig[eglsCHKPT] = 1;
1393 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1398 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1400 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1402 gmx_bool bIfRandomize;
1403 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1404 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1405 if (constr && bIfRandomize)
1407 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1408 state, fr->bMolPBC, graph, f,
1409 &top->idef, tmp_vir,
1410 cr, nrnb, wcycle, upd, constr,
1411 TRUE, bCalcVir, vetanew);
1416 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1418 gmx_iterate_init(&iterate, TRUE);
1419 /* for iterations, we save these vectors, as we will be redoing the calculations */
1420 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1423 bFirstIterate = TRUE;
1424 while (bFirstIterate || iterate.bIterationActive)
1426 /* We now restore these vectors to redo the calculation with improved extended variables */
1427 if (iterate.bIterationActive)
1429 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1432 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1433 so scroll down for that logic */
1435 /* ######### START SECOND UPDATE STEP ################# */
1436 /* Box is changed in update() when we do pressure coupling,
1437 * but we should still use the old box for energy corrections and when
1438 * writing it to the energy file, so it matches the trajectory files for
1439 * the same timestep above. Make a copy in a separate array.
1441 copy_mat(state->box, lastbox);
1446 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1448 wallcycle_start(wcycle, ewcUPDATE);
1449 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1452 if (iterate.bIterationActive)
1460 /* we use a new value of scalevir to converge the iterations faster */
1461 scalevir = tracevir/trace(shake_vir);
1463 msmul(shake_vir, scalevir, shake_vir);
1464 m_add(force_vir, shake_vir, total_vir);
1465 clear_mat(shake_vir);
1467 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1468 /* We can only do Berendsen coupling after we have summed
1469 * the kinetic energy or virial. Since the happens
1470 * in global_state after update, we should only do it at
1471 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1476 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1477 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1482 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1484 /* velocity half-step update */
1485 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1486 bUpdateDoLR, fr->f_twin, fcd,
1487 ekind, M, upd, FALSE, etrtVELOCITY2,
1488 cr, nrnb, constr, &top->idef);
1491 /* Above, initialize just copies ekinh into ekin,
1492 * it doesn't copy position (for VV),
1493 * and entire integrator for MD.
1496 if (ir->eI == eiVVAK)
1498 copy_rvecn(state->x, cbuf, 0, state->natoms);
1500 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1502 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1503 bUpdateDoLR, fr->f_twin, fcd,
1504 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1505 wallcycle_stop(wcycle, ewcUPDATE);
1507 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1508 fr->bMolPBC, graph, f,
1509 &top->idef, shake_vir,
1510 cr, nrnb, wcycle, upd, constr,
1511 FALSE, bCalcVir, state->veta);
1513 if (ir->eI == eiVVAK)
1515 /* erase F_EKIN and F_TEMP here? */
1516 /* just compute the kinetic energy at the half step to perform a trotter step */
1517 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1518 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1519 constr, NULL, FALSE, lastbox,
1520 top_global, &bSumEkinhOld,
1521 cglo_flags | CGLO_TEMPERATURE
1523 wallcycle_start(wcycle, ewcUPDATE);
1524 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1525 /* now we know the scaling, we can compute the positions again again */
1526 copy_rvecn(cbuf, state->x, 0, state->natoms);
1528 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1530 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1531 bUpdateDoLR, fr->f_twin, fcd,
1532 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1533 wallcycle_stop(wcycle, ewcUPDATE);
1535 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1536 /* are the small terms in the shake_vir here due
1537 * to numerical errors, or are they important
1538 * physically? I'm thinking they are just errors, but not completely sure.
1539 * For now, will call without actually constraining, constr=NULL*/
1540 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1541 state, fr->bMolPBC, graph, f,
1542 &top->idef, tmp_vir,
1543 cr, nrnb, wcycle, upd, NULL,
1549 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1552 if (fr->bSepDVDL && fplog && do_log)
1554 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1558 /* this factor or 2 correction is necessary
1559 because half of the constraint force is removed
1560 in the vv step, so we have to double it. See
1561 the Redmine issue #1255. It is not yet clear
1562 if the factor of 2 is exact, or just a very
1563 good approximation, and this will be
1564 investigated. The next step is to see if this
1565 can be done adding a dhdl contribution from the
1566 rattle step, but this is somewhat more
1567 complicated with the current code. Will be
1568 investigated, hopefully for 4.6.3. However,
1569 this current solution is much better than
1570 having it completely wrong.
1572 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1576 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1581 /* Need to unshift here */
1582 unshift_self(graph, state->box, state->x);
1587 wallcycle_start(wcycle, ewcVSITECONSTR);
1590 shift_self(graph, state->box, state->x);
1592 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1593 top->idef.iparams, top->idef.il,
1594 fr->ePBC, fr->bMolPBC, cr, state->box);
1598 unshift_self(graph, state->box, state->x);
1600 wallcycle_stop(wcycle, ewcVSITECONSTR);
1603 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1604 /* With Leap-Frog we can skip compute_globals at
1605 * non-communication steps, but we need to calculate
1606 * the kinetic energy one step before communication.
1608 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1610 if (ir->nstlist == -1 && bFirstIterate)
1612 gs.sig[eglsNABNSB] = nlh.nabnsb;
1614 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1615 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1617 bFirstIterate ? &gs : NULL,
1618 (step_rel % gs.nstms == 0) &&
1619 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1621 top_global, &bSumEkinhOld,
1623 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1624 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1625 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1626 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1627 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1628 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1631 if (ir->nstlist == -1 && bFirstIterate)
1633 nlh.nabnsb = gs.set[eglsNABNSB];
1634 gs.set[eglsNABNSB] = 0;
1637 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1638 /* ############# END CALC EKIN AND PRESSURE ################# */
1640 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1641 the virial that should probably be addressed eventually. state->veta has better properies,
1642 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1643 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1645 if (iterate.bIterationActive &&
1646 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1647 trace(shake_vir), &tracevir))
1651 bFirstIterate = FALSE;
1654 if (!bVV || bRerunMD)
1656 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1657 sum_dhdl(enerd, state->lambda, ir->fepvals);
1659 update_box(fplog, step, ir, mdatoms, state, f,
1660 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1662 /* ################# END UPDATE STEP 2 ################# */
1663 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1665 /* The coordinates (x) were unshifted in update */
1668 /* We will not sum ekinh_old,
1669 * so signal that we still have to do it.
1671 bSumEkinhOld = TRUE;
1674 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1676 /* use the directly determined last velocity, not actually the averaged half steps */
1677 if (bTrotter && ir->eI == eiVV)
1679 enerd->term[F_EKIN] = last_ekin;
1681 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1685 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1689 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1691 /* ######### END PREPARING EDR OUTPUT ########### */
1696 gmx_bool do_dr, do_or;
1698 if (fplog && do_log && bDoExpanded)
1700 /* only needed if doing expanded ensemble */
1701 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1702 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1704 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1708 upd_mdebin(mdebin, bDoDHDL, TRUE,
1709 t, mdatoms->tmass, enerd, state,
1710 ir->fepvals, ir->expandedvals, lastbox,
1711 shake_vir, force_vir, total_vir, pres,
1712 ekind, mu_tot, constr);
1716 upd_mdebin_step(mdebin);
1719 do_dr = do_per_step(step, ir->nstdisreout);
1720 do_or = do_per_step(step, ir->nstorireout);
1722 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1724 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1726 if (ir->ePull != epullNO)
1728 pull_print_output(ir->pull, step, t);
1731 if (do_per_step(step, ir->nstlog))
1733 if (fflush(fplog) != 0)
1735 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1741 /* Have to do this part _after_ outputting the logfile and the edr file */
1742 /* Gets written into the state at the beginning of next loop*/
1743 state->fep_state = lamnew;
1745 /* Print the remaining wall clock time for the run */
1746 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1750 fprintf(stderr, "\n");
1752 print_time(stderr, walltime_accounting, step, ir, cr);
1755 /* Ion/water position swapping.
1756 * Not done in last step since trajectory writing happens before this call
1757 * in the MD loop and exchanges would be lost anyway. */
1758 bNeedRepartition = FALSE;
1759 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1760 do_per_step(step, ir->swap->nstswap))
1762 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1763 bRerunMD ? rerun_fr.x : state->x,
1764 bRerunMD ? rerun_fr.box : state->box,
1765 top_global, MASTER(cr) && bVerbose, bRerunMD);
1767 if (bNeedRepartition && DOMAINDECOMP(cr))
1769 dd_collect_state(cr->dd, state, state_global);
1773 /* Replica exchange */
1775 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1776 do_per_step(step, repl_ex_nst))
1778 bExchanged = replica_exchange(fplog, cr, repl_ex,
1779 state_global, enerd,
1783 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1785 dd_partition_system(fplog, step, cr, TRUE, 1,
1786 state_global, top_global, ir,
1787 state, &f, mdatoms, top, fr,
1788 vsite, shellfc, constr,
1789 nrnb, wcycle, FALSE);
1794 bStartingFromCpt = FALSE;
1796 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1797 /* With all integrators, except VV, we need to retain the pressure
1798 * at the current step for coupling at the next step.
1800 if ((state->flags & (1<<estPRES_PREV)) &&
1802 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1804 /* Store the pressure in t_state for pressure coupling
1805 * at the next MD step.
1807 copy_mat(pres, state->pres_prev);
1810 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1812 if ( (membed != NULL) && (!bLastStep) )
1814 rescale_membed(step_rel, membed, state_global->x);
1821 /* read next frame from input trajectory */
1822 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1827 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1831 if (!bRerunMD || !rerun_fr.bStep)
1833 /* increase the MD step number */
1838 cycles = wallcycle_stop(wcycle, ewcSTEP);
1839 if (DOMAINDECOMP(cr) && wcycle)
1841 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1844 if (bPMETuneRunning || bPMETuneTry)
1846 /* PME grid + cut-off optimization with GPUs or PME nodes */
1848 /* Count the total cycles over the last steps */
1849 cycles_pmes += cycles;
1851 /* We can only switch cut-off at NS steps */
1852 if (step % ir->nstlist == 0)
1854 /* PME grid + cut-off optimization with GPUs or PME nodes */
1857 if (DDMASTER(cr->dd))
1859 /* PME node load is too high, start tuning */
1860 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1862 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1864 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1866 bPMETuneTry = FALSE;
1869 if (bPMETuneRunning)
1871 /* init_step might not be a multiple of nstlist,
1872 * but the first cycle is always skipped anyhow.
1875 pme_load_balance(pme_loadbal, cr,
1876 (bVerbose && MASTER(cr)) ? stderr : NULL,
1878 ir, state, cycles_pmes,
1879 fr->ic, fr->nbv, &fr->pmedata,
1882 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1883 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1884 fr->ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1885 fr->rlist = fr->ic->rlist;
1886 fr->rlistlong = fr->ic->rlistlong;
1887 fr->rcoulomb = fr->ic->rcoulomb;
1888 fr->rvdw = fr->ic->rvdw;
1890 if (ir->eDispCorr != edispcNO)
1892 calc_enervirdiff(NULL, ir->eDispCorr, fr);
1899 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1900 gs.set[eglsRESETCOUNTERS] != 0)
1902 /* Reset all the counters related to performance over the run */
1903 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1904 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1905 wcycle_set_reset_counters(wcycle, -1);
1906 if (!(cr->duty & DUTY_PME))
1908 /* Tell our PME node to reset its counters */
1909 gmx_pme_send_resetcounters(cr, step);
1911 /* Correct max_hours for the elapsed time */
1912 max_hours -= elapsed_time/(60.0*60.0);
1913 bResetCountersHalfMaxH = FALSE;
1914 gs.set[eglsRESETCOUNTERS] = 0;
1918 /* End of main MD loop */
1921 /* Stop measuring walltime */
1922 walltime_accounting_end(walltime_accounting);
1924 if (bRerunMD && MASTER(cr))
1929 if (!(cr->duty & DUTY_PME))
1931 /* Tell the PME only node to finish */
1932 gmx_pme_send_finish(cr);
1937 if (ir->nstcalcenergy > 0 && !bRerunMD)
1939 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1940 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1947 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1949 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)));
1950 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1953 if (pme_loadbal != NULL)
1955 pme_loadbal_done(pme_loadbal, cr, fplog,
1956 fr->nbv != NULL && fr->nbv->bUseGPU);
1959 if (shellfc && fplog)
1961 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1962 (nconverged*100.0)/step_rel);
1963 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1967 if (repl_ex_nst > 0 && MASTER(cr))
1969 print_replica_exchange_statistics(fplog, repl_ex);
1972 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);