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 if ((Flags & MD_TUNEPME) &&
463 EEL_PME(fr->eeltype) &&
464 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
465 !bRerunMD && !EVDW_PME(fr->vdwtype))
467 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
469 if (cr->duty & DUTY_PME)
471 /* Start tuning right away, as we can't measure the load */
472 bPMETuneRunning = TRUE;
476 /* Separate PME nodes, we can measure the PP/PME load balance */
481 if (!ir->bContinuation && !bRerunMD)
483 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
485 /* Set the velocities of frozen particles to zero */
486 for (i = 0; i < mdatoms->homenr; i++)
488 for (m = 0; m < DIM; m++)
490 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
500 /* Constrain the initial coordinates and velocities */
501 do_constrain_first(fplog, constr, ir, mdatoms, state,
506 /* Construct the virtual sites for the initial configuration */
507 construct_vsites(vsite, state->x, ir->delta_t, NULL,
508 top->idef.iparams, top->idef.il,
509 fr->ePBC, fr->bMolPBC, cr, state->box);
515 /* set free energy calculation frequency as the minimum
516 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
517 nstfep = ir->fepvals->nstdhdl;
520 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
524 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
527 /* I'm assuming we need global communication the first time! MRS */
528 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
529 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
530 | (bVV ? CGLO_PRESSURE : 0)
531 | (bVV ? CGLO_CONSTRAINT : 0)
532 | (bRerunMD ? CGLO_RERUNMD : 0)
533 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
535 bSumEkinhOld = FALSE;
536 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
537 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
538 constr, NULL, FALSE, state->box,
539 top_global, &bSumEkinhOld, cglo_flags);
540 if (ir->eI == eiVVAK)
542 /* a second call to get the half step temperature initialized as well */
543 /* we do the same call as above, but turn the pressure off -- internally to
544 compute_globals, this is recognized as a velocity verlet half-step
545 kinetic energy calculation. This minimized excess variables, but
546 perhaps loses some logic?*/
548 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
549 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
550 constr, NULL, FALSE, state->box,
551 top_global, &bSumEkinhOld,
552 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
555 /* Calculate the initial half step temperature, and save the ekinh_old */
556 if (!(Flags & MD_STARTFROMCPT))
558 for (i = 0; (i < ir->opts.ngtc); i++)
560 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
565 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
566 and there is no previous step */
569 /* if using an iterative algorithm, we need to create a working directory for the state. */
572 bufstate = init_bufstate(state);
575 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
576 temperature control */
577 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
581 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
584 "RMS relative constraint deviation after constraining: %.2e\n",
585 constr_rmsd(constr, FALSE));
587 if (EI_STATE_VELOCITY(ir->eI))
589 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
593 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
594 " input trajectory '%s'\n\n",
595 *(top_global->name), opt2fn("-rerun", nfile, fnm));
598 fprintf(stderr, "Calculated time to finish depends on nsteps from "
599 "run input file,\nwhich may not correspond to the time "
600 "needed to process input trajectory.\n\n");
606 fprintf(stderr, "starting mdrun '%s'\n",
607 *(top_global->name));
610 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
614 sprintf(tbuf, "%s", "infinite");
616 if (ir->init_step > 0)
618 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
619 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
620 gmx_step_str(ir->init_step, sbuf2),
621 ir->init_step*ir->delta_t);
625 fprintf(stderr, "%s steps, %s ps.\n",
626 gmx_step_str(ir->nsteps, sbuf), tbuf);
629 fprintf(fplog, "\n");
632 walltime_accounting_start(walltime_accounting);
633 wallcycle_start(wcycle, ewcRUN);
634 print_start(fplog, cr, walltime_accounting, "mdrun");
636 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
638 chkpt_ret = fcCheckPointParallel( cr->nodeid,
642 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
647 /***********************************************************
651 ************************************************************/
653 /* if rerunMD then read coordinates and velocities from input trajectory */
656 if (getenv("GMX_FORCE_UPDATE"))
664 bNotLastFrame = read_first_frame(oenv, &status,
665 opt2fn("-rerun", nfile, fnm),
666 &rerun_fr, TRX_NEED_X | TRX_READ_V);
667 if (rerun_fr.natoms != top_global->natoms)
670 "Number of atoms in trajectory (%d) does not match the "
671 "run input file (%d)\n",
672 rerun_fr.natoms, top_global->natoms);
674 if (ir->ePBC != epbcNONE)
678 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);
680 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
682 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
689 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
692 if (ir->ePBC != epbcNONE)
694 /* Set the shift vectors.
695 * Necessary here when have a static box different from the tpr box.
697 calc_shifts(rerun_fr.box, fr->shift_vec);
701 /* loop over MD steps or if rerunMD to end of input trajectory */
703 /* Skip the first Nose-Hoover integration when we get the state from tpx */
704 bStateFromTPX = !bStateFromCP;
705 bInitStep = bFirstStep && (bStateFromTPX || bVV);
706 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
708 bSumEkinhOld = FALSE;
710 bNeedRepartition = FALSE;
712 init_global_signals(&gs, cr, ir, repl_ex_nst);
714 step = ir->init_step;
717 if (ir->nstlist == -1)
719 init_nlistheuristics(&nlh, bGStatEveryStep, step);
722 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
724 /* check how many steps are left in other sims */
725 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
729 /* and stop now if we should */
730 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
731 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
732 while (!bLastStep || (bRerunMD && bNotLastFrame))
735 wallcycle_start(wcycle, ewcSTEP);
741 step = rerun_fr.step;
742 step_rel = step - ir->init_step;
755 bLastStep = (step_rel == ir->nsteps);
756 t = t0 + step*ir->delta_t;
759 if (ir->efep != efepNO || ir->bSimTemp)
761 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
762 requiring different logic. */
764 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
765 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
766 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
767 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
768 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
773 update_annealing_target_temp(&(ir->opts), t);
778 if (!DOMAINDECOMP(cr) || MASTER(cr))
780 for (i = 0; i < state_global->natoms; i++)
782 copy_rvec(rerun_fr.x[i], state_global->x[i]);
786 for (i = 0; i < state_global->natoms; i++)
788 copy_rvec(rerun_fr.v[i], state_global->v[i]);
793 for (i = 0; i < state_global->natoms; i++)
795 clear_rvec(state_global->v[i]);
799 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
800 " Ekin, temperature and pressure are incorrect,\n"
801 " the virial will be incorrect when constraints are present.\n"
803 bRerunWarnNoV = FALSE;
807 copy_mat(rerun_fr.box, state_global->box);
808 copy_mat(state_global->box, state->box);
810 if (vsite && (Flags & MD_RERUN_VSITE))
812 if (DOMAINDECOMP(cr))
814 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
818 /* Following is necessary because the graph may get out of sync
819 * with the coordinates if we only have every N'th coordinate set
821 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
822 shift_self(graph, state->box, state->x);
824 construct_vsites(vsite, state->x, ir->delta_t, state->v,
825 top->idef.iparams, top->idef.il,
826 fr->ePBC, fr->bMolPBC, cr, state->box);
829 unshift_self(graph, state->box, state->x);
834 /* Stop Center of Mass motion */
835 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
839 /* for rerun MD always do Neighbour Searching */
840 bNS = (bFirstStep || ir->nstlist != 0);
845 /* Determine whether or not to do Neighbour Searching and LR */
846 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
848 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
849 (ir->nstlist == -1 && nlh.nabnsb > 0));
851 if (bNS && ir->nstlist == -1)
853 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
857 /* check whether we should stop because another simulation has
861 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
862 (multisim_nsteps != ir->nsteps) )
869 "Stopping simulation %d because another one has finished\n",
873 gs.sig[eglsCHKPT] = 1;
878 /* < 0 means stop at next step, > 0 means stop at next NS step */
879 if ( (gs.set[eglsSTOPCOND] < 0) ||
880 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
885 /* Determine whether or not to update the Born radii if doing GB */
886 bBornRadii = bFirstStep;
887 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
892 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
893 do_verbose = bVerbose &&
894 (step % stepout == 0 || bFirstStep || bLastStep);
896 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
904 bMasterState = FALSE;
905 /* Correct the new box if it is too skewed */
906 if (DYNAMIC_BOX(*ir))
908 if (correct_box(fplog, step, state->box, graph))
913 if (DOMAINDECOMP(cr) && bMasterState)
915 dd_collect_state(cr->dd, state, state_global);
919 if (DOMAINDECOMP(cr))
921 /* Repartition the domain decomposition */
922 wallcycle_start(wcycle, ewcDOMDEC);
923 dd_partition_system(fplog, step, cr,
924 bMasterState, nstglobalcomm,
925 state_global, top_global, ir,
926 state, &f, mdatoms, top, fr,
927 vsite, shellfc, constr,
929 do_verbose && !bPMETuneRunning);
930 wallcycle_stop(wcycle, ewcDOMDEC);
931 /* If using an iterative integrator, reallocate space to match the decomposition */
935 if (MASTER(cr) && do_log)
937 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
940 if (ir->efep != efepNO)
942 update_mdatoms(mdatoms, state->lambda[efptMASS]);
945 if ((bRerunMD && rerun_fr.bV) || bExchanged)
948 /* We need the kinetic energy at minus the half step for determining
949 * the full step kinetic energy and possibly for T-coupling.*/
950 /* This may not be quite working correctly yet . . . . */
951 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
952 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
953 constr, NULL, FALSE, state->box,
954 top_global, &bSumEkinhOld,
955 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
957 clear_mat(force_vir);
959 /* We write a checkpoint at this MD step when:
960 * either at an NS step when we signalled through gs,
961 * or at the last step (but not when we do not want confout),
962 * but never at the first step or with rerun.
964 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
965 (bLastStep && (Flags & MD_CONFOUT))) &&
966 step > ir->init_step && !bRerunMD);
969 gs.set[eglsCHKPT] = 0;
972 /* Determine the energy and pressure:
973 * at nstcalcenergy steps and at energy output steps (set below).
975 if (EI_VV(ir->eI) && (!bInitStep))
977 /* for vv, the first half of the integration actually corresponds
978 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
979 but the virial needs to be calculated on both the current step and the 'next' step. Future
980 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
982 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
983 bCalcVir = bCalcEner ||
984 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
988 bCalcEner = do_per_step(step, ir->nstcalcenergy);
989 bCalcVir = bCalcEner ||
990 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
993 /* Do we need global communication ? */
994 bGStat = (bCalcVir || bCalcEner || bStopCM ||
995 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
996 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
998 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1000 if (do_ene || do_log)
1007 /* these CGLO_ options remain the same throughout the iteration */
1008 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1009 (bGStat ? CGLO_GSTAT : 0)
1012 force_flags = (GMX_FORCE_STATECHANGED |
1013 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1014 GMX_FORCE_ALLFORCES |
1016 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1017 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1018 (bDoFEP ? GMX_FORCE_DHDL : 0)
1023 if (do_per_step(step, ir->nstcalclr))
1025 force_flags |= GMX_FORCE_DO_LR;
1031 /* Now is the time to relax the shells */
1032 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1033 ir, bNS, force_flags,
1036 state, f, force_vir, mdatoms,
1037 nrnb, wcycle, graph, groups,
1038 shellfc, fr, bBornRadii, t, mu_tot,
1040 mdoutf_get_fp_field(outf));
1050 /* The coordinates (x) are shifted (to get whole molecules)
1052 * This is parallellized as well, and does communication too.
1053 * Check comments in sim_util.c
1055 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1056 state->box, state->x, &state->hist,
1057 f, force_vir, mdatoms, enerd, fcd,
1058 state->lambda, graph,
1059 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1060 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1063 if (bVV && !bStartingFromCpt && !bRerunMD)
1064 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1066 if (ir->eI == eiVV && bInitStep)
1068 /* if using velocity verlet with full time step Ekin,
1069 * take the first half step only to compute the
1070 * virial for the first step. From there,
1071 * revert back to the initial coordinates
1072 * so that the input is actually the initial step.
1074 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1078 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1079 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1082 /* If we are using twin-range interactions where the long-range component
1083 * is only evaluated every nstcalclr>1 steps, we should do a special update
1084 * step to combine the long-range forces on these steps.
1085 * For nstcalclr=1 this is not done, since the forces would have been added
1086 * directly to the short-range forces already.
1088 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1090 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1091 f, bUpdateDoLR, fr->f_twin, fcd,
1092 ekind, M, upd, bInitStep, etrtVELOCITY1,
1093 cr, nrnb, constr, &top->idef);
1095 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1097 gmx_iterate_init(&iterate, TRUE);
1099 /* for iterations, we save these vectors, as we will be self-consistently iterating
1102 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1104 /* save the state */
1105 if (iterate.bIterationActive)
1107 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1110 bFirstIterate = TRUE;
1111 while (bFirstIterate || iterate.bIterationActive)
1113 if (iterate.bIterationActive)
1115 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1116 if (bFirstIterate && bTrotter)
1118 /* The first time through, we need a decent first estimate
1119 of veta(t+dt) to compute the constraints. Do
1120 this by computing the box volume part of the
1121 trotter integration at this time. Nothing else
1122 should be changed by this routine here. If
1123 !(first time), we start with the previous value
1126 veta_save = state->veta;
1127 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1128 vetanew = state->veta;
1129 state->veta = veta_save;
1134 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1136 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1137 state, fr->bMolPBC, graph, f,
1138 &top->idef, shake_vir,
1139 cr, nrnb, wcycle, upd, constr,
1140 TRUE, bCalcVir, vetanew);
1144 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1150 /* Need to unshift here if a do_force has been
1151 called in the previous step */
1152 unshift_self(graph, state->box, state->x);
1155 /* if VV, compute the pressure and constraints */
1156 /* For VV2, we strictly only need this if using pressure
1157 * control, but we really would like to have accurate pressures
1159 * Think about ways around this in the future?
1160 * For now, keep this choice in comments.
1162 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1163 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1165 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1166 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1168 bSumEkinhOld = TRUE;
1170 /* for vv, the first half of the integration actually corresponds to the previous step.
1171 So we need information from the last step in the first half of the integration */
1172 if (bGStat || do_per_step(step-1, nstglobalcomm))
1174 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1175 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1176 constr, NULL, FALSE, state->box,
1177 top_global, &bSumEkinhOld,
1180 | (bTemp ? CGLO_TEMPERATURE : 0)
1181 | (bPres ? CGLO_PRESSURE : 0)
1182 | (bPres ? CGLO_CONSTRAINT : 0)
1183 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1184 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1187 /* explanation of above:
1188 a) We compute Ekin at the full time step
1189 if 1) we are using the AveVel Ekin, and it's not the
1190 initial step, or 2) if we are using AveEkin, but need the full
1191 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1192 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1193 EkinAveVel because it's needed for the pressure */
1195 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1200 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1201 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1208 /* We need the kinetic energy at minus the half step for determining
1209 * the full step kinetic energy and possibly for T-coupling.*/
1210 /* This may not be quite working correctly yet . . . . */
1211 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1212 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1213 constr, NULL, FALSE, state->box,
1214 top_global, &bSumEkinhOld,
1215 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1220 if (iterate.bIterationActive &&
1221 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1222 state->veta, &vetanew))
1226 bFirstIterate = FALSE;
1229 if (bTrotter && !bInitStep)
1231 copy_mat(shake_vir, state->svir_prev);
1232 copy_mat(force_vir, state->fvir_prev);
1233 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1235 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1236 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1237 enerd->term[F_EKIN] = trace(ekind->ekin);
1240 /* if it's the initial step, we performed this first step just to get the constraint virial */
1241 if (bInitStep && ir->eI == eiVV)
1243 copy_rvecn(cbuf, state->v, 0, state->natoms);
1247 /* MRS -- now done iterating -- compute the conserved quantity */
1250 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1253 last_ekin = enerd->term[F_EKIN];
1255 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1257 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1259 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1262 sum_dhdl(enerd, state->lambda, ir->fepvals);
1266 /* ######## END FIRST UPDATE STEP ############## */
1267 /* ######## If doing VV, we now have v(dt) ###### */
1270 /* perform extended ensemble sampling in lambda - we don't
1271 actually move to the new state before outputting
1272 statistics, but if performing simulated tempering, we
1273 do update the velocities and the tau_t. */
1275 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1276 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1277 copy_df_history(&state_global->dfhist, &state->dfhist);
1280 /* Now we have the energies and forces corresponding to the
1281 * coordinates at time t. We must output all of this before
1284 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1285 ir, state, state_global, top_global, fr,
1286 outf, mdebin, ekind, f, f_global,
1288 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1291 /* kludge -- virial is lost with restart for NPT control. Must restart */
1292 if (bStartingFromCpt && bVV)
1294 copy_mat(state->svir_prev, shake_vir);
1295 copy_mat(state->fvir_prev, force_vir);
1298 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1300 /* Check whether everything is still allright */
1301 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1302 #ifdef GMX_THREAD_MPI
1307 /* this is just make gs.sig compatible with the hack
1308 of sending signals around by MPI_Reduce with together with
1310 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1312 gs.sig[eglsSTOPCOND] = 1;
1314 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1316 gs.sig[eglsSTOPCOND] = -1;
1318 /* < 0 means stop at next step, > 0 means stop at next NS step */
1322 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1323 gmx_get_signal_name(),
1324 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1328 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1329 gmx_get_signal_name(),
1330 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1332 handled_stop_condition = (int)gmx_get_stop_condition();
1334 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1335 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1336 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1338 /* Signal to terminate the run */
1339 gs.sig[eglsSTOPCOND] = 1;
1342 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1344 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1347 if (bResetCountersHalfMaxH && MASTER(cr) &&
1348 elapsed_time > max_hours*60.0*60.0*0.495)
1350 gs.sig[eglsRESETCOUNTERS] = 1;
1353 if (ir->nstlist == -1 && !bRerunMD)
1355 /* When bGStatEveryStep=FALSE, global_stat is only called
1356 * when we check the atom displacements, not at NS steps.
1357 * This means that also the bonded interaction count check is not
1358 * performed immediately after NS. Therefore a few MD steps could
1359 * be performed with missing interactions.
1360 * But wrong energies are never written to file,
1361 * since energies are only written after global_stat
1364 if (step >= nlh.step_nscheck)
1366 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1367 nlh.scale_tot, state->x);
1371 /* This is not necessarily true,
1372 * but step_nscheck is determined quite conservatively.
1378 /* In parallel we only have to check for checkpointing in steps
1379 * where we do global communication,
1380 * otherwise the other nodes don't know.
1382 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1385 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1386 gs.set[eglsCHKPT] == 0)
1388 gs.sig[eglsCHKPT] = 1;
1391 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1396 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1398 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1400 gmx_bool bIfRandomize;
1401 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1402 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1403 if (constr && bIfRandomize)
1405 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1406 state, fr->bMolPBC, graph, f,
1407 &top->idef, tmp_vir,
1408 cr, nrnb, wcycle, upd, constr,
1409 TRUE, bCalcVir, vetanew);
1414 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1416 gmx_iterate_init(&iterate, TRUE);
1417 /* for iterations, we save these vectors, as we will be redoing the calculations */
1418 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1421 bFirstIterate = TRUE;
1422 while (bFirstIterate || iterate.bIterationActive)
1424 /* We now restore these vectors to redo the calculation with improved extended variables */
1425 if (iterate.bIterationActive)
1427 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1430 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1431 so scroll down for that logic */
1433 /* ######### START SECOND UPDATE STEP ################# */
1434 /* Box is changed in update() when we do pressure coupling,
1435 * but we should still use the old box for energy corrections and when
1436 * writing it to the energy file, so it matches the trajectory files for
1437 * the same timestep above. Make a copy in a separate array.
1439 copy_mat(state->box, lastbox);
1444 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1446 wallcycle_start(wcycle, ewcUPDATE);
1447 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1450 if (iterate.bIterationActive)
1458 /* we use a new value of scalevir to converge the iterations faster */
1459 scalevir = tracevir/trace(shake_vir);
1461 msmul(shake_vir, scalevir, shake_vir);
1462 m_add(force_vir, shake_vir, total_vir);
1463 clear_mat(shake_vir);
1465 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1466 /* We can only do Berendsen coupling after we have summed
1467 * the kinetic energy or virial. Since the happens
1468 * in global_state after update, we should only do it at
1469 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1474 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1475 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1480 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1482 /* velocity half-step update */
1483 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1484 bUpdateDoLR, fr->f_twin, fcd,
1485 ekind, M, upd, FALSE, etrtVELOCITY2,
1486 cr, nrnb, constr, &top->idef);
1489 /* Above, initialize just copies ekinh into ekin,
1490 * it doesn't copy position (for VV),
1491 * and entire integrator for MD.
1494 if (ir->eI == eiVVAK)
1496 copy_rvecn(state->x, cbuf, 0, state->natoms);
1498 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1500 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1501 bUpdateDoLR, fr->f_twin, fcd,
1502 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1503 wallcycle_stop(wcycle, ewcUPDATE);
1505 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1506 fr->bMolPBC, graph, f,
1507 &top->idef, shake_vir,
1508 cr, nrnb, wcycle, upd, constr,
1509 FALSE, bCalcVir, state->veta);
1511 if (ir->eI == eiVVAK)
1513 /* erase F_EKIN and F_TEMP here? */
1514 /* just compute the kinetic energy at the half step to perform a trotter step */
1515 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1516 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1517 constr, NULL, FALSE, lastbox,
1518 top_global, &bSumEkinhOld,
1519 cglo_flags | CGLO_TEMPERATURE
1521 wallcycle_start(wcycle, ewcUPDATE);
1522 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1523 /* now we know the scaling, we can compute the positions again again */
1524 copy_rvecn(cbuf, state->x, 0, state->natoms);
1526 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1528 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1529 bUpdateDoLR, fr->f_twin, fcd,
1530 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1531 wallcycle_stop(wcycle, ewcUPDATE);
1533 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1534 /* are the small terms in the shake_vir here due
1535 * to numerical errors, or are they important
1536 * physically? I'm thinking they are just errors, but not completely sure.
1537 * For now, will call without actually constraining, constr=NULL*/
1538 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1539 state, fr->bMolPBC, graph, f,
1540 &top->idef, tmp_vir,
1541 cr, nrnb, wcycle, upd, NULL,
1547 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1550 if (fr->bSepDVDL && fplog && do_log)
1552 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1556 /* this factor or 2 correction is necessary
1557 because half of the constraint force is removed
1558 in the vv step, so we have to double it. See
1559 the Redmine issue #1255. It is not yet clear
1560 if the factor of 2 is exact, or just a very
1561 good approximation, and this will be
1562 investigated. The next step is to see if this
1563 can be done adding a dhdl contribution from the
1564 rattle step, but this is somewhat more
1565 complicated with the current code. Will be
1566 investigated, hopefully for 4.6.3. However,
1567 this current solution is much better than
1568 having it completely wrong.
1570 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1574 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1579 /* Need to unshift here */
1580 unshift_self(graph, state->box, state->x);
1585 wallcycle_start(wcycle, ewcVSITECONSTR);
1588 shift_self(graph, state->box, state->x);
1590 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1591 top->idef.iparams, top->idef.il,
1592 fr->ePBC, fr->bMolPBC, cr, state->box);
1596 unshift_self(graph, state->box, state->x);
1598 wallcycle_stop(wcycle, ewcVSITECONSTR);
1601 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1602 /* With Leap-Frog we can skip compute_globals at
1603 * non-communication steps, but we need to calculate
1604 * the kinetic energy one step before communication.
1606 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1608 if (ir->nstlist == -1 && bFirstIterate)
1610 gs.sig[eglsNABNSB] = nlh.nabnsb;
1612 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1613 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1615 bFirstIterate ? &gs : NULL,
1616 (step_rel % gs.nstms == 0) &&
1617 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1619 top_global, &bSumEkinhOld,
1621 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1622 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1623 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1624 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1625 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1626 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1629 if (ir->nstlist == -1 && bFirstIterate)
1631 nlh.nabnsb = gs.set[eglsNABNSB];
1632 gs.set[eglsNABNSB] = 0;
1635 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1636 /* ############# END CALC EKIN AND PRESSURE ################# */
1638 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1639 the virial that should probably be addressed eventually. state->veta has better properies,
1640 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1641 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1643 if (iterate.bIterationActive &&
1644 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1645 trace(shake_vir), &tracevir))
1649 bFirstIterate = FALSE;
1652 if (!bVV || bRerunMD)
1654 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1655 sum_dhdl(enerd, state->lambda, ir->fepvals);
1657 update_box(fplog, step, ir, mdatoms, state, f,
1658 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1660 /* ################# END UPDATE STEP 2 ################# */
1661 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1663 /* The coordinates (x) were unshifted in update */
1666 /* We will not sum ekinh_old,
1667 * so signal that we still have to do it.
1669 bSumEkinhOld = TRUE;
1672 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1674 /* use the directly determined last velocity, not actually the averaged half steps */
1675 if (bTrotter && ir->eI == eiVV)
1677 enerd->term[F_EKIN] = last_ekin;
1679 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1683 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1687 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1689 /* ######### END PREPARING EDR OUTPUT ########### */
1694 gmx_bool do_dr, do_or;
1696 if (fplog && do_log && bDoExpanded)
1698 /* only needed if doing expanded ensemble */
1699 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1700 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1702 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1706 upd_mdebin(mdebin, bDoDHDL, TRUE,
1707 t, mdatoms->tmass, enerd, state,
1708 ir->fepvals, ir->expandedvals, lastbox,
1709 shake_vir, force_vir, total_vir, pres,
1710 ekind, mu_tot, constr);
1714 upd_mdebin_step(mdebin);
1717 do_dr = do_per_step(step, ir->nstdisreout);
1718 do_or = do_per_step(step, ir->nstorireout);
1720 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1722 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1724 if (ir->ePull != epullNO)
1726 pull_print_output(ir->pull, step, t);
1729 if (do_per_step(step, ir->nstlog))
1731 if (fflush(fplog) != 0)
1733 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1739 /* Have to do this part _after_ outputting the logfile and the edr file */
1740 /* Gets written into the state at the beginning of next loop*/
1741 state->fep_state = lamnew;
1743 /* Print the remaining wall clock time for the run */
1744 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1748 fprintf(stderr, "\n");
1750 print_time(stderr, walltime_accounting, step, ir, cr);
1753 /* Ion/water position swapping.
1754 * Not done in last step since trajectory writing happens before this call
1755 * in the MD loop and exchanges would be lost anyway. */
1756 bNeedRepartition = FALSE;
1757 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1758 do_per_step(step, ir->swap->nstswap))
1760 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1761 bRerunMD ? rerun_fr.x : state->x,
1762 bRerunMD ? rerun_fr.box : state->box,
1763 top_global, MASTER(cr) && bVerbose, bRerunMD);
1765 if (bNeedRepartition && DOMAINDECOMP(cr))
1767 dd_collect_state(cr->dd, state, state_global);
1771 /* Replica exchange */
1773 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1774 do_per_step(step, repl_ex_nst))
1776 bExchanged = replica_exchange(fplog, cr, repl_ex,
1777 state_global, enerd,
1781 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1783 dd_partition_system(fplog, step, cr, TRUE, 1,
1784 state_global, top_global, ir,
1785 state, &f, mdatoms, top, fr,
1786 vsite, shellfc, constr,
1787 nrnb, wcycle, FALSE);
1792 bStartingFromCpt = FALSE;
1794 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1795 /* With all integrators, except VV, we need to retain the pressure
1796 * at the current step for coupling at the next step.
1798 if ((state->flags & (1<<estPRES_PREV)) &&
1800 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1802 /* Store the pressure in t_state for pressure coupling
1803 * at the next MD step.
1805 copy_mat(pres, state->pres_prev);
1808 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1810 if ( (membed != NULL) && (!bLastStep) )
1812 rescale_membed(step_rel, membed, state_global->x);
1819 /* read next frame from input trajectory */
1820 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1825 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1829 if (!bRerunMD || !rerun_fr.bStep)
1831 /* increase the MD step number */
1836 cycles = wallcycle_stop(wcycle, ewcSTEP);
1837 if (DOMAINDECOMP(cr) && wcycle)
1839 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1842 if (bPMETuneRunning || bPMETuneTry)
1844 /* PME grid + cut-off optimization with GPUs or PME nodes */
1846 /* Count the total cycles over the last steps */
1847 cycles_pmes += cycles;
1849 /* We can only switch cut-off at NS steps */
1850 if (step % ir->nstlist == 0)
1852 /* PME grid + cut-off optimization with GPUs or PME nodes */
1855 if (DDMASTER(cr->dd))
1857 /* PME node load is too high, start tuning */
1858 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1860 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1862 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1864 bPMETuneTry = FALSE;
1867 if (bPMETuneRunning)
1869 /* init_step might not be a multiple of nstlist,
1870 * but the first cycle is always skipped anyhow.
1873 pme_load_balance(pme_loadbal, cr,
1874 (bVerbose && MASTER(cr)) ? stderr : NULL,
1876 ir, state, cycles_pmes,
1877 fr->ic, fr->nbv, &fr->pmedata,
1880 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1881 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1882 fr->rlist = fr->ic->rlist;
1883 fr->rlistlong = fr->ic->rlistlong;
1884 fr->rcoulomb = fr->ic->rcoulomb;
1885 fr->rvdw = fr->ic->rvdw;
1891 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1892 gs.set[eglsRESETCOUNTERS] != 0)
1894 /* Reset all the counters related to performance over the run */
1895 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1896 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1897 wcycle_set_reset_counters(wcycle, -1);
1898 if (!(cr->duty & DUTY_PME))
1900 /* Tell our PME node to reset its counters */
1901 gmx_pme_send_resetcounters(cr, step);
1903 /* Correct max_hours for the elapsed time */
1904 max_hours -= elapsed_time/(60.0*60.0);
1905 bResetCountersHalfMaxH = FALSE;
1906 gs.set[eglsRESETCOUNTERS] = 0;
1910 /* End of main MD loop */
1913 /* Stop measuring walltime */
1914 walltime_accounting_end(walltime_accounting);
1916 if (bRerunMD && MASTER(cr))
1921 if (!(cr->duty & DUTY_PME))
1923 /* Tell the PME only node to finish */
1924 gmx_pme_send_finish(cr);
1929 if (ir->nstcalcenergy > 0 && !bRerunMD)
1931 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1932 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1939 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1941 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)));
1942 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1945 if (pme_loadbal != NULL)
1947 pme_loadbal_done(pme_loadbal, cr, fplog,
1948 fr->nbv != NULL && fr->nbv->bUseGPU);
1951 if (shellfc && fplog)
1953 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1954 (nconverged*100.0)/step_rel);
1955 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1959 if (repl_ex_nst > 0 && MASTER(cr))
1961 print_replica_exchange_statistics(fplog, repl_ex);
1964 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);