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,2015, 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.
45 #include "thread_mpi/threads.h"
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/domdec/domdec_network.h"
49 #include "gromacs/ewald/pme-load-balancing.h"
50 #include "gromacs/ewald/pme.h"
51 #include "gromacs/fileio/filenm.h"
52 #include "gromacs/fileio/mdoutf.h"
53 #include "gromacs/fileio/trajectory_writing.h"
54 #include "gromacs/fileio/trx.h"
55 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/imd/imd.h"
57 #include "gromacs/legacyheaders/constr.h"
58 #include "gromacs/legacyheaders/ebin.h"
59 #include "gromacs/legacyheaders/force.h"
60 #include "gromacs/legacyheaders/md_logging.h"
61 #include "gromacs/legacyheaders/md_support.h"
62 #include "gromacs/legacyheaders/mdatoms.h"
63 #include "gromacs/legacyheaders/mdebin.h"
64 #include "gromacs/legacyheaders/mdrun.h"
65 #include "gromacs/legacyheaders/network.h"
66 #include "gromacs/legacyheaders/nrnb.h"
67 #include "gromacs/legacyheaders/ns.h"
68 #include "gromacs/legacyheaders/shellfc.h"
69 #include "gromacs/legacyheaders/sighandler.h"
70 #include "gromacs/legacyheaders/sim_util.h"
71 #include "gromacs/legacyheaders/tgroup.h"
72 #include "gromacs/legacyheaders/typedefs.h"
73 #include "gromacs/legacyheaders/update.h"
74 #include "gromacs/legacyheaders/vcm.h"
75 #include "gromacs/legacyheaders/vsite.h"
76 #include "gromacs/legacyheaders/types/commrec.h"
77 #include "gromacs/legacyheaders/types/constr.h"
78 #include "gromacs/legacyheaders/types/enums.h"
79 #include "gromacs/legacyheaders/types/fcdata.h"
80 #include "gromacs/legacyheaders/types/force_flags.h"
81 #include "gromacs/legacyheaders/types/forcerec.h"
82 #include "gromacs/legacyheaders/types/group.h"
83 #include "gromacs/legacyheaders/types/inputrec.h"
84 #include "gromacs/legacyheaders/types/interaction_const.h"
85 #include "gromacs/legacyheaders/types/mdatom.h"
86 #include "gromacs/legacyheaders/types/membedt.h"
87 #include "gromacs/legacyheaders/types/nrnb.h"
88 #include "gromacs/legacyheaders/types/oenv.h"
89 #include "gromacs/legacyheaders/types/shellfc.h"
90 #include "gromacs/legacyheaders/types/state.h"
91 #include "gromacs/listed-forces/manage-threading.h"
92 #include "gromacs/math/utilities.h"
93 #include "gromacs/math/vec.h"
94 #include "gromacs/math/vectypes.h"
95 #include "gromacs/mdlib/compute_io.h"
96 #include "gromacs/mdlib/mdrun_signalling.h"
97 #include "gromacs/mdlib/nb_verlet.h"
98 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
99 #include "gromacs/pbcutil/mshift.h"
100 #include "gromacs/pbcutil/pbc.h"
101 #include "gromacs/pulling/pull.h"
102 #include "gromacs/swap/swapcoords.h"
103 #include "gromacs/timing/wallcycle.h"
104 #include "gromacs/timing/walltime_accounting.h"
105 #include "gromacs/topology/atoms.h"
106 #include "gromacs/topology/idef.h"
107 #include "gromacs/topology/mtop_util.h"
108 #include "gromacs/topology/topology.h"
109 #include "gromacs/utility/basedefinitions.h"
110 #include "gromacs/utility/cstringutil.h"
111 #include "gromacs/utility/fatalerror.h"
112 #include "gromacs/utility/real.h"
113 #include "gromacs/utility/smalloc.h"
120 #include "corewrap.h"
123 static void reset_all_counters(FILE *fplog, t_commrec *cr,
125 gmx_int64_t *step_rel, t_inputrec *ir,
126 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
127 gmx_walltime_accounting_t walltime_accounting,
128 struct nonbonded_verlet_t *nbv)
130 char sbuf[STEPSTRSIZE];
132 /* Reset all the counters related to performance over the run */
133 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
134 gmx_step_str(step, sbuf));
138 nbnxn_gpu_reset_timings(nbv);
141 wallcycle_stop(wcycle, ewcRUN);
142 wallcycle_reset_all(wcycle);
143 if (DOMAINDECOMP(cr))
145 reset_dd_statistics_counters(cr->dd);
148 ir->init_step += *step_rel;
149 ir->nsteps -= *step_rel;
151 wallcycle_start(wcycle, ewcRUN);
152 walltime_accounting_start(walltime_accounting);
153 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
156 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
157 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
159 gmx_vsite_t *vsite, gmx_constr_t constr,
160 int stepout, t_inputrec *ir,
161 gmx_mtop_t *top_global,
163 t_state *state_global,
165 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
166 gmx_edsam_t ed, t_forcerec *fr,
167 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
168 real cpt_period, real max_hours,
171 gmx_walltime_accounting_t walltime_accounting)
173 gmx_mdoutf_t outf = NULL;
174 gmx_int64_t step, step_rel;
176 double t, t0, lam0[efptNR];
177 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
178 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
179 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
180 bBornRadii, bStartingFromCpt;
181 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
182 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
183 bForceUpdate = FALSE, bCPT;
184 gmx_bool bMasterState;
185 int force_flags, cglo_flags;
186 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
193 gmx_repl_ex_t repl_ex = NULL;
196 t_mdebin *mdebin = NULL;
197 t_state *state = NULL;
198 rvec *f_global = NULL;
199 gmx_enerdata_t *enerd;
201 gmx_global_stat_t gstat;
202 gmx_update_t upd = NULL;
203 t_graph *graph = NULL;
205 gmx_groups_t *groups;
206 gmx_ekindata_t *ekind;
207 gmx_shellfc_t shellfc;
208 int count, nconverged = 0;
210 gmx_bool bConverged = TRUE, bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
211 gmx_bool bResetCountersHalfMaxH = FALSE;
212 gmx_bool bVV, bTemp, bPres, bTrotter;
213 gmx_bool bUpdateDoLR;
222 real saved_conserved_quantity = 0;
226 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
227 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
228 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
229 simulation stops. If equal to zero, don't
230 communicate any more between multisims.*/
231 /* PME load balancing data for GPU kernels */
232 pme_load_balancing_t *pme_loadbal;
233 gmx_bool bPMETune = FALSE;
234 gmx_bool bPMETunePrinting = FALSE;
237 gmx_bool bIMDstep = FALSE;
240 /* Temporary addition for FAHCORE checkpointing */
244 /* Check for special mdrun options */
245 bRerunMD = (Flags & MD_RERUN);
246 if (Flags & MD_RESETCOUNTERSHALFWAY)
250 /* Signal to reset the counters half the simulation steps. */
251 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
253 /* Signal to reset the counters halfway the simulation time. */
254 bResetCountersHalfMaxH = (max_hours > 0);
257 /* md-vv uses averaged full step velocities for T-control
258 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
259 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
261 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
265 /* Since we don't know if the frames read are related in any way,
266 * rebuild the neighborlist at every step.
269 ir->nstcalcenergy = 1;
273 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
275 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
276 bGStatEveryStep = (nstglobalcomm == 1);
280 ir->nstxout_compressed = 0;
282 groups = &top_global->groups;
285 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
286 &(state_global->fep_state), lam0,
287 nrnb, top_global, &upd,
288 nfile, fnm, &outf, &mdebin,
289 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags, wcycle);
291 clear_mat(total_vir);
293 /* Energy terms and groups */
295 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
297 if (DOMAINDECOMP(cr))
303 snew(f, top_global->natoms);
306 /* Kinetic energy data */
308 init_ekindata(fplog, top_global, &(ir->opts), ekind);
309 /* Copy the cos acceleration to the groups struct */
310 ekind->cosacc.cos_accel = ir->cos_accel;
312 gstat = global_stat_init(ir);
315 /* Check for polarizable models and flexible constraints */
316 shellfc = init_shell_flexcon(fplog,
317 top_global, n_flexible_constraints(constr),
318 (ir->bContinuation ||
319 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
320 NULL : state_global->x);
321 if (shellfc && ir->nstcalcenergy != 1)
323 gmx_fatal(FARGS, "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is not supported in combinations with shell particles.\nPlease make a new tpr file.", ir->nstcalcenergy);
325 if (shellfc && DOMAINDECOMP(cr))
327 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
329 if (shellfc && ir->eI == eiNM)
331 /* Currently shells don't work with Normal Modes */
332 gmx_fatal(FARGS, "Normal Mode analysis is not supported with shells.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
335 if (vsite && ir->eI == eiNM)
337 /* Currently virtual sites don't work with Normal Modes */
338 gmx_fatal(FARGS, "Normal Mode analysis is not supported with virtual sites.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
341 if (bRerunMD && fr->cutoff_scheme == ecutsVERLET && ir->opts.ngener > 1 && usingGpu(fr->nbv))
343 gmx_fatal(FARGS, "The Verlet scheme on GPUs does not support energy groups, so your rerun should probably use a .tpr file without energy groups, or mdrun -nb auto");
348 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
349 set_deform_reference_box(upd,
350 deform_init_init_step_tpx,
351 deform_init_box_tpx);
352 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
356 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
357 if ((io > 2000) && MASTER(cr))
360 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
365 if (DOMAINDECOMP(cr))
367 top = dd_init_local_top(top_global);
370 dd_init_local_state(cr->dd, state_global, state);
372 if (DDMASTER(cr->dd) && ir->nstfout)
374 snew(f_global, state_global->natoms);
379 top = gmx_mtop_generate_local_top(top_global, ir);
381 forcerec_set_excl_load(fr, top);
383 state = serial_init_local_state(state_global);
386 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
390 set_vsite_top(vsite, top, mdatoms, cr);
393 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
395 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
400 make_local_shells(cr, mdatoms, shellfc);
403 setup_bonded_threading(fr, &top->idef);
406 /* Set up interactive MD (IMD) */
407 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
408 nfile, fnm, oenv, imdport, Flags);
410 if (DOMAINDECOMP(cr))
412 /* Distribute the charge groups over the nodes from the master node */
413 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
414 state_global, top_global, ir,
415 state, &f, mdatoms, top, fr,
416 vsite, shellfc, constr,
421 update_mdatoms(mdatoms, state->lambda[efptMASS]);
423 if (opt2bSet("-cpi", nfile, fnm))
425 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
429 bStateFromCP = FALSE;
434 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
441 /* Update mdebin with energy history if appending to output files */
442 if (Flags & MD_APPENDFILES)
444 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
448 /* We might have read an energy history from checkpoint,
449 * free the allocated memory and reset the counts.
451 done_energyhistory(&state_global->enerhist);
452 init_energyhistory(&state_global->enerhist);
455 /* Set the initial energy history in state by updating once */
456 update_energyhistory(&state_global->enerhist, mdebin);
459 /* Initialize constraints */
460 if (constr && !DOMAINDECOMP(cr))
462 set_constraints(constr, top, ir, mdatoms, cr);
465 if (repl_ex_nst > 0 && MASTER(cr))
467 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
468 repl_ex_nst, repl_ex_nex, repl_ex_seed);
471 /* PME tuning is only supported with PME for Coulomb. Is is not supported
472 * with only LJ PME, or for reruns.
474 bPMETune = ((Flags & MD_TUNEPME) && EEL_PME(fr->eeltype) && !bRerunMD &&
475 !(Flags & MD_REPRODUCIBLE));
478 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata,
479 use_GPU(fr->nbv), !(cr->duty & DUTY_PME),
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 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
572 temperature control */
573 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
577 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
580 "RMS relative constraint deviation after constraining: %.2e\n",
581 constr_rmsd(constr, FALSE));
583 if (EI_STATE_VELOCITY(ir->eI))
585 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
589 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
590 " input trajectory '%s'\n\n",
591 *(top_global->name), opt2fn("-rerun", nfile, fnm));
594 fprintf(stderr, "Calculated time to finish depends on nsteps from "
595 "run input file,\nwhich may not correspond to the time "
596 "needed to process input trajectory.\n\n");
602 fprintf(stderr, "starting mdrun '%s'\n",
603 *(top_global->name));
606 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
610 sprintf(tbuf, "%s", "infinite");
612 if (ir->init_step > 0)
614 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
615 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
616 gmx_step_str(ir->init_step, sbuf2),
617 ir->init_step*ir->delta_t);
621 fprintf(stderr, "%s steps, %s ps.\n",
622 gmx_step_str(ir->nsteps, sbuf), tbuf);
625 fprintf(fplog, "\n");
628 walltime_accounting_start(walltime_accounting);
629 wallcycle_start(wcycle, ewcRUN);
630 print_start(fplog, cr, walltime_accounting, "mdrun");
632 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
634 chkpt_ret = fcCheckPointParallel( cr->nodeid,
638 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
643 /***********************************************************
647 ************************************************************/
649 /* if rerunMD then read coordinates and velocities from input trajectory */
652 if (getenv("GMX_FORCE_UPDATE"))
660 bNotLastFrame = read_first_frame(oenv, &status,
661 opt2fn("-rerun", nfile, fnm),
662 &rerun_fr, TRX_NEED_X | TRX_READ_V);
663 if (rerun_fr.natoms != top_global->natoms)
666 "Number of atoms in trajectory (%d) does not match the "
667 "run input file (%d)\n",
668 rerun_fr.natoms, top_global->natoms);
670 if (ir->ePBC != epbcNONE)
674 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);
676 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
678 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
685 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
688 if (ir->ePBC != epbcNONE)
690 /* Set the shift vectors.
691 * Necessary here when have a static box different from the tpr box.
693 calc_shifts(rerun_fr.box, fr->shift_vec);
697 /* loop over MD steps or if rerunMD to end of input trajectory */
699 /* Skip the first Nose-Hoover integration when we get the state from tpx */
700 bStateFromTPX = !bStateFromCP;
701 bInitStep = bFirstStep && (bStateFromTPX || bVV);
702 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
703 bSumEkinhOld = FALSE;
705 bNeedRepartition = FALSE;
707 init_global_signals(&gs, cr, ir, repl_ex_nst);
709 step = ir->init_step;
712 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
714 /* check how many steps are left in other sims */
715 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
719 /* and stop now if we should */
720 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
721 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
722 while (!bLastStep || (bRerunMD && bNotLastFrame))
725 /* Determine if this is a neighbor search step */
726 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
728 if (bPMETune && bNStList)
730 /* PME grid + cut-off optimization with GPUs or PME nodes */
731 pme_loadbal_do(pme_loadbal, cr,
732 (bVerbose && MASTER(cr)) ? stderr : NULL,
734 ir, fr, state, wcycle,
739 wallcycle_start(wcycle, ewcSTEP);
745 step = rerun_fr.step;
746 step_rel = step - ir->init_step;
759 bLastStep = (step_rel == ir->nsteps);
760 t = t0 + step*ir->delta_t;
763 if (ir->efep != efepNO || ir->bSimTemp)
765 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
766 requiring different logic. */
768 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
769 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
770 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
771 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
772 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
775 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
776 do_per_step(step, repl_ex_nst));
780 update_annealing_target_temp(&(ir->opts), t);
785 if (!DOMAINDECOMP(cr) || MASTER(cr))
787 for (i = 0; i < state_global->natoms; i++)
789 copy_rvec(rerun_fr.x[i], state_global->x[i]);
793 for (i = 0; i < state_global->natoms; i++)
795 copy_rvec(rerun_fr.v[i], state_global->v[i]);
800 for (i = 0; i < state_global->natoms; i++)
802 clear_rvec(state_global->v[i]);
806 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
807 " Ekin, temperature and pressure are incorrect,\n"
808 " the virial will be incorrect when constraints are present.\n"
810 bRerunWarnNoV = FALSE;
814 copy_mat(rerun_fr.box, state_global->box);
815 copy_mat(state_global->box, state->box);
817 if (vsite && (Flags & MD_RERUN_VSITE))
819 if (DOMAINDECOMP(cr))
821 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
825 /* Following is necessary because the graph may get out of sync
826 * with the coordinates if we only have every N'th coordinate set
828 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
829 shift_self(graph, state->box, state->x);
831 construct_vsites(vsite, state->x, ir->delta_t, state->v,
832 top->idef.iparams, top->idef.il,
833 fr->ePBC, fr->bMolPBC, cr, state->box);
836 unshift_self(graph, state->box, state->x);
841 /* Stop Center of Mass motion */
842 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
846 /* for rerun MD always do Neighbour Searching */
847 bNS = (bFirstStep || ir->nstlist != 0);
852 bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP);
855 /* check whether we should stop because another simulation has
859 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
860 (multisim_nsteps != ir->nsteps) )
867 "Stopping simulation %d because another one has finished\n",
871 gs.sig[eglsCHKPT] = 1;
876 /* < 0 means stop at next step, > 0 means stop at next NS step */
877 if ( (gs.set[eglsSTOPCOND] < 0) ||
878 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
883 /* Determine whether or not to update the Born radii if doing GB */
884 bBornRadii = bFirstStep;
885 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
890 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
891 do_verbose = bVerbose &&
892 (step % stepout == 0 || bFirstStep || bLastStep);
894 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
902 bMasterState = FALSE;
903 /* Correct the new box if it is too skewed */
904 if (DYNAMIC_BOX(*ir))
906 if (correct_box(fplog, step, state->box, graph))
911 if (DOMAINDECOMP(cr) && bMasterState)
913 dd_collect_state(cr->dd, state, state_global);
917 if (DOMAINDECOMP(cr))
919 /* Repartition the domain decomposition */
920 dd_partition_system(fplog, step, cr,
921 bMasterState, nstglobalcomm,
922 state_global, top_global, ir,
923 state, &f, mdatoms, top, fr,
924 vsite, shellfc, constr,
926 do_verbose && !bPMETunePrinting);
930 if (MASTER(cr) && do_log)
932 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
935 if (ir->efep != efepNO)
937 update_mdatoms(mdatoms, state->lambda[efptMASS]);
940 if ((bRerunMD && rerun_fr.bV) || bExchanged)
943 /* We need the kinetic energy at minus the half step for determining
944 * the full step kinetic energy and possibly for T-coupling.*/
945 /* This may not be quite working correctly yet . . . . */
946 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
947 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
948 constr, NULL, FALSE, state->box,
949 top_global, &bSumEkinhOld,
950 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
952 clear_mat(force_vir);
954 /* We write a checkpoint at this MD step when:
955 * either at an NS step when we signalled through gs,
956 * or at the last step (but not when we do not want confout),
957 * but never at the first step or with rerun.
959 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
960 (bLastStep && (Flags & MD_CONFOUT))) &&
961 step > ir->init_step && !bRerunMD);
964 gs.set[eglsCHKPT] = 0;
967 /* Determine the energy and pressure:
968 * at nstcalcenergy steps and at energy output steps (set below).
970 if (EI_VV(ir->eI) && (!bInitStep))
972 /* for vv, the first half of the integration actually corresponds
973 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
974 but the virial needs to be calculated on both the current step and the 'next' step. Future
975 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
977 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
978 bCalcVir = bCalcEner ||
979 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
983 bCalcEner = do_per_step(step, ir->nstcalcenergy);
984 bCalcVir = bCalcEner ||
985 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
988 /* Do we need global communication ? */
989 bGStat = (bCalcVir || bCalcEner || bStopCM ||
990 do_per_step(step, nstglobalcomm) ||
991 (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)));
993 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
995 if (do_ene || do_log || bDoReplEx)
1002 /* these CGLO_ options remain the same throughout the iteration */
1003 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1004 (bGStat ? CGLO_GSTAT : 0)
1007 force_flags = (GMX_FORCE_STATECHANGED |
1008 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1009 GMX_FORCE_ALLFORCES |
1011 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1012 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1013 (bDoFEP ? GMX_FORCE_DHDL : 0)
1018 if (do_per_step(step, ir->nstcalclr))
1020 force_flags |= GMX_FORCE_DO_LR;
1026 /* Now is the time to relax the shells */
1027 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1028 ir, bNS, force_flags,
1031 state, f, force_vir, mdatoms,
1032 nrnb, wcycle, graph, groups,
1033 shellfc, fr, bBornRadii, t, mu_tot,
1035 mdoutf_get_fp_field(outf));
1045 /* The coordinates (x) are shifted (to get whole molecules)
1047 * This is parallellized as well, and does communication too.
1048 * Check comments in sim_util.c
1050 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1051 state->box, state->x, &state->hist,
1052 f, force_vir, mdatoms, enerd, fcd,
1053 state->lambda, graph,
1054 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1055 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1058 if (bVV && !bStartingFromCpt && !bRerunMD)
1059 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1063 wallcycle_start(wcycle, ewcUPDATE);
1064 if (ir->eI == eiVV && bInitStep)
1066 /* if using velocity verlet with full time step Ekin,
1067 * take the first half step only to compute the
1068 * virial for the first step. From there,
1069 * revert back to the initial coordinates
1070 * so that the input is actually the initial step.
1072 snew(vbuf, state->natoms);
1073 copy_rvecn(state->v, vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1077 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1078 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1081 /* If we are using twin-range interactions where the long-range component
1082 * is only evaluated every nstcalclr>1 steps, we should do a special update
1083 * step to combine the long-range forces on these steps.
1084 * For nstcalclr=1 this is not done, since the forces would have been added
1085 * directly to the short-range forces already.
1087 * TODO Remove various aspects of VV+twin-range in master
1088 * branch, because VV integrators did not ever support
1089 * twin-range multiple time stepping with constraints.
1091 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1093 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1094 f, bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1095 ekind, M, upd, bInitStep, etrtVELOCITY1,
1096 cr, nrnb, constr, &top->idef);
1098 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1100 wallcycle_stop(wcycle, ewcUPDATE);
1101 update_constraints(fplog, step, NULL, ir, mdatoms,
1102 state, fr->bMolPBC, graph, f,
1103 &top->idef, shake_vir,
1104 cr, nrnb, wcycle, upd, constr,
1106 wallcycle_start(wcycle, ewcUPDATE);
1107 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1109 /* Correct the virial for multiple time stepping */
1110 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1115 /* Need to unshift here if a do_force has been
1116 called in the previous step */
1117 unshift_self(graph, state->box, state->x);
1119 /* if VV, compute the pressure and constraints */
1120 /* For VV2, we strictly only need this if using pressure
1121 * control, but we really would like to have accurate pressures
1123 * Think about ways around this in the future?
1124 * For now, keep this choice in comments.
1126 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1127 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1129 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1130 if (bCalcEner && ir->eI == eiVVAK)
1132 bSumEkinhOld = TRUE;
1134 /* for vv, the first half of the integration actually corresponds to the previous step.
1135 So we need information from the last step in the first half of the integration */
1136 if (bGStat || do_per_step(step-1, nstglobalcomm))
1138 wallcycle_stop(wcycle, ewcUPDATE);
1139 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1140 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1141 constr, NULL, FALSE, state->box,
1142 top_global, &bSumEkinhOld,
1145 | (bTemp ? CGLO_TEMPERATURE : 0)
1146 | (bPres ? CGLO_PRESSURE : 0)
1147 | (bPres ? CGLO_CONSTRAINT : 0)
1150 /* explanation of above:
1151 a) We compute Ekin at the full time step
1152 if 1) we are using the AveVel Ekin, and it's not the
1153 initial step, or 2) if we are using AveEkin, but need the full
1154 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1155 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1156 EkinAveVel because it's needed for the pressure */
1157 wallcycle_start(wcycle, ewcUPDATE);
1159 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1164 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1165 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1171 wallcycle_stop(wcycle, ewcUPDATE);
1172 /* We need the kinetic energy at minus the half step for determining
1173 * the full step kinetic energy and possibly for T-coupling.*/
1174 /* This may not be quite working correctly yet . . . . */
1175 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1176 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1177 constr, NULL, FALSE, state->box,
1178 top_global, &bSumEkinhOld,
1179 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1180 wallcycle_start(wcycle, ewcUPDATE);
1184 if (bTrotter && !bInitStep)
1186 copy_mat(shake_vir, state->svir_prev);
1187 copy_mat(force_vir, state->fvir_prev);
1188 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1190 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1191 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1192 enerd->term[F_EKIN] = trace(ekind->ekin);
1195 /* if it's the initial step, we performed this first step just to get the constraint virial */
1196 if (ir->eI == eiVV && bInitStep)
1198 copy_rvecn(vbuf, state->v, 0, state->natoms);
1201 wallcycle_stop(wcycle, ewcUPDATE);
1204 /* compute the conserved quantity */
1207 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1210 last_ekin = enerd->term[F_EKIN];
1212 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1214 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1216 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1217 if (ir->efep != efepNO && !bRerunMD)
1219 sum_dhdl(enerd, state->lambda, ir->fepvals);
1223 /* ######## END FIRST UPDATE STEP ############## */
1224 /* ######## If doing VV, we now have v(dt) ###### */
1227 /* perform extended ensemble sampling in lambda - we don't
1228 actually move to the new state before outputting
1229 statistics, but if performing simulated tempering, we
1230 do update the velocities and the tau_t. */
1232 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1233 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1234 copy_df_history(&state_global->dfhist, &state->dfhist);
1237 /* Now we have the energies and forces corresponding to the
1238 * coordinates at time t. We must output all of this before
1241 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1242 ir, state, state_global, top_global, fr,
1243 outf, mdebin, ekind, f, f_global,
1245 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1247 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1248 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1250 /* kludge -- virial is lost with restart for NPT control. Must restart */
1251 if (bStartingFromCpt && bVV)
1253 copy_mat(state->svir_prev, shake_vir);
1254 copy_mat(state->fvir_prev, force_vir);
1257 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1259 /* Check whether everything is still allright */
1260 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1261 #ifdef GMX_THREAD_MPI
1266 /* this is just make gs.sig compatible with the hack
1267 of sending signals around by MPI_Reduce with together with
1269 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1271 gs.sig[eglsSTOPCOND] = 1;
1273 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1275 gs.sig[eglsSTOPCOND] = -1;
1277 /* < 0 means stop at next step, > 0 means stop at next NS step */
1281 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1282 gmx_get_signal_name(),
1283 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1287 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1288 gmx_get_signal_name(),
1289 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1291 handled_stop_condition = (int)gmx_get_stop_condition();
1293 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1294 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1295 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1297 /* Signal to terminate the run */
1298 gs.sig[eglsSTOPCOND] = 1;
1301 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1303 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1306 if (bResetCountersHalfMaxH && MASTER(cr) &&
1307 elapsed_time > max_hours*60.0*60.0*0.495)
1309 gs.sig[eglsRESETCOUNTERS] = 1;
1312 /* In parallel we only have to check for checkpointing in steps
1313 * where we do global communication,
1314 * otherwise the other nodes don't know.
1316 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1319 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1320 gs.set[eglsCHKPT] == 0)
1322 gs.sig[eglsCHKPT] = 1;
1325 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1330 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1332 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1334 gmx_bool bIfRandomize;
1335 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1336 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1337 if (constr && bIfRandomize)
1339 update_constraints(fplog, step, NULL, ir, mdatoms,
1340 state, fr->bMolPBC, graph, f,
1341 &top->idef, tmp_vir,
1342 cr, nrnb, wcycle, upd, constr,
1347 /* ######### START SECOND UPDATE STEP ################# */
1348 /* Box is changed in update() when we do pressure coupling,
1349 * but we should still use the old box for energy corrections and when
1350 * writing it to the energy file, so it matches the trajectory files for
1351 * the same timestep above. Make a copy in a separate array.
1353 copy_mat(state->box, lastbox);
1357 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1359 wallcycle_start(wcycle, ewcUPDATE);
1360 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1363 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1364 /* We can only do Berendsen coupling after we have summed
1365 * the kinetic energy or virial. Since the happens
1366 * in global_state after update, we should only do it at
1367 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1372 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1373 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1378 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1380 /* velocity half-step update */
1381 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1382 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1383 ekind, M, upd, FALSE, etrtVELOCITY2,
1384 cr, nrnb, constr, &top->idef);
1387 /* Above, initialize just copies ekinh into ekin,
1388 * it doesn't copy position (for VV),
1389 * and entire integrator for MD.
1392 if (ir->eI == eiVVAK)
1394 /* We probably only need md->homenr, not state->natoms */
1395 if (state->natoms > cbuf_nalloc)
1397 cbuf_nalloc = state->natoms;
1398 srenew(cbuf, cbuf_nalloc);
1400 copy_rvecn(state->x, cbuf, 0, state->natoms);
1402 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1404 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1405 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1406 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1407 wallcycle_stop(wcycle, ewcUPDATE);
1409 update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1410 fr->bMolPBC, graph, f,
1411 &top->idef, shake_vir,
1412 cr, nrnb, wcycle, upd, constr,
1415 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1417 /* Correct the virial for multiple time stepping */
1418 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1421 if (ir->eI == eiVVAK)
1423 /* erase F_EKIN and F_TEMP here? */
1424 /* just compute the kinetic energy at the half step to perform a trotter step */
1425 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1426 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1427 constr, NULL, FALSE, lastbox,
1428 top_global, &bSumEkinhOld,
1429 cglo_flags | CGLO_TEMPERATURE
1431 wallcycle_start(wcycle, ewcUPDATE);
1432 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1433 /* now we know the scaling, we can compute the positions again again */
1434 copy_rvecn(cbuf, state->x, 0, state->natoms);
1436 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1438 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1439 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1440 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1441 wallcycle_stop(wcycle, ewcUPDATE);
1443 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1444 /* are the small terms in the shake_vir here due
1445 * to numerical errors, or are they important
1446 * physically? I'm thinking they are just errors, but not completely sure.
1447 * For now, will call without actually constraining, constr=NULL*/
1448 update_constraints(fplog, step, NULL, ir, mdatoms,
1449 state, fr->bMolPBC, graph, f,
1450 &top->idef, tmp_vir,
1451 cr, nrnb, wcycle, upd, NULL,
1456 /* this factor or 2 correction is necessary
1457 because half of the constraint force is removed
1458 in the vv step, so we have to double it. See
1459 the Redmine issue #1255. It is not yet clear
1460 if the factor of 2 is exact, or just a very
1461 good approximation, and this will be
1462 investigated. The next step is to see if this
1463 can be done adding a dhdl contribution from the
1464 rattle step, but this is somewhat more
1465 complicated with the current code. Will be
1466 investigated, hopefully for 4.6.3. However,
1467 this current solution is much better than
1468 having it completely wrong.
1470 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1474 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1479 /* Need to unshift here */
1480 unshift_self(graph, state->box, state->x);
1485 wallcycle_start(wcycle, ewcVSITECONSTR);
1488 shift_self(graph, state->box, state->x);
1490 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1491 top->idef.iparams, top->idef.il,
1492 fr->ePBC, fr->bMolPBC, cr, state->box);
1496 unshift_self(graph, state->box, state->x);
1498 wallcycle_stop(wcycle, ewcVSITECONSTR);
1501 /* ############## IF NOT VV, Calculate globals HERE ############ */
1502 /* With Leap-Frog we can skip compute_globals at
1503 * non-communication steps, but we need to calculate
1504 * the kinetic energy one step before communication.
1506 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1508 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1509 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1511 (step_rel % gs.nstms == 0) &&
1512 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1514 top_global, &bSumEkinhOld,
1516 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1517 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1518 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1519 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1524 /* ############# END CALC EKIN AND PRESSURE ################# */
1526 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1527 the virial that should probably be addressed eventually. state->veta has better properies,
1528 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1529 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1531 if (ir->efep != efepNO && (!bVV || bRerunMD))
1533 /* Sum up the foreign energy and dhdl terms for md and sd.
1534 Currently done every step so that dhdl is correct in the .edr */
1535 sum_dhdl(enerd, state->lambda, ir->fepvals);
1537 update_box(fplog, step, ir, mdatoms, state, f,
1538 pcoupl_mu, nrnb, upd);
1540 /* ################# END UPDATE STEP 2 ################# */
1541 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1543 /* The coordinates (x) were unshifted in update */
1546 /* We will not sum ekinh_old,
1547 * so signal that we still have to do it.
1549 bSumEkinhOld = TRUE;
1552 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1554 /* use the directly determined last velocity, not actually the averaged half steps */
1555 if (bTrotter && ir->eI == eiVV)
1557 enerd->term[F_EKIN] = last_ekin;
1559 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1563 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1567 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1569 /* ######### END PREPARING EDR OUTPUT ########### */
1574 gmx_bool do_dr, do_or;
1576 if (fplog && do_log && bDoExpanded)
1578 /* only needed if doing expanded ensemble */
1579 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1580 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1582 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1586 upd_mdebin(mdebin, bDoDHDL, TRUE,
1587 t, mdatoms->tmass, enerd, state,
1588 ir->fepvals, ir->expandedvals, lastbox,
1589 shake_vir, force_vir, total_vir, pres,
1590 ekind, mu_tot, constr);
1594 upd_mdebin_step(mdebin);
1597 do_dr = do_per_step(step, ir->nstdisreout);
1598 do_or = do_per_step(step, ir->nstorireout);
1600 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1602 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1606 pull_print_output(ir->pull, step, t);
1609 if (do_per_step(step, ir->nstlog))
1611 if (fflush(fplog) != 0)
1613 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1619 /* Have to do this part _after_ outputting the logfile and the edr file */
1620 /* Gets written into the state at the beginning of next loop*/
1621 state->fep_state = lamnew;
1623 /* Print the remaining wall clock time for the run */
1624 if (MULTIMASTER(cr) &&
1625 (do_verbose || gmx_got_usr_signal()) &&
1630 fprintf(stderr, "\n");
1632 print_time(stderr, walltime_accounting, step, ir, cr);
1635 /* Ion/water position swapping.
1636 * Not done in last step since trajectory writing happens before this call
1637 * in the MD loop and exchanges would be lost anyway. */
1638 bNeedRepartition = FALSE;
1639 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1640 do_per_step(step, ir->swap->nstswap))
1642 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1643 bRerunMD ? rerun_fr.x : state->x,
1644 bRerunMD ? rerun_fr.box : state->box,
1645 top_global, MASTER(cr) && bVerbose, bRerunMD);
1647 if (bNeedRepartition && DOMAINDECOMP(cr))
1649 dd_collect_state(cr->dd, state, state_global);
1653 /* Replica exchange */
1657 bExchanged = replica_exchange(fplog, cr, repl_ex,
1658 state_global, enerd,
1662 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1664 dd_partition_system(fplog, step, cr, TRUE, 1,
1665 state_global, top_global, ir,
1666 state, &f, mdatoms, top, fr,
1667 vsite, shellfc, constr,
1668 nrnb, wcycle, FALSE);
1673 bStartingFromCpt = FALSE;
1675 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1676 /* With all integrators, except VV, we need to retain the pressure
1677 * at the current step for coupling at the next step.
1679 if ((state->flags & (1<<estPRES_PREV)) &&
1681 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1683 /* Store the pressure in t_state for pressure coupling
1684 * at the next MD step.
1686 copy_mat(pres, state->pres_prev);
1689 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1691 if ( (membed != NULL) && (!bLastStep) )
1693 rescale_membed(step_rel, membed, state_global->x);
1700 /* read next frame from input trajectory */
1701 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1706 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1710 cycles = wallcycle_stop(wcycle, ewcSTEP);
1711 if (DOMAINDECOMP(cr) && wcycle)
1713 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1716 if (!bRerunMD || !rerun_fr.bStep)
1718 /* increase the MD step number */
1723 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1724 gs.set[eglsRESETCOUNTERS] != 0)
1726 /* Reset all the counters related to performance over the run */
1727 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1728 use_GPU(fr->nbv) ? fr->nbv : NULL);
1729 wcycle_set_reset_counters(wcycle, -1);
1730 if (!(cr->duty & DUTY_PME))
1732 /* Tell our PME node to reset its counters */
1733 gmx_pme_send_resetcounters(cr, step);
1735 /* Correct max_hours for the elapsed time */
1736 max_hours -= elapsed_time/(60.0*60.0);
1737 bResetCountersHalfMaxH = FALSE;
1738 gs.set[eglsRESETCOUNTERS] = 0;
1741 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1742 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1745 /* End of main MD loop */
1748 /* Closing TNG files can include compressing data. Therefore it is good to do that
1749 * before stopping the time measurements. */
1750 mdoutf_tng_close(outf);
1752 /* Stop measuring walltime */
1753 walltime_accounting_end(walltime_accounting);
1755 if (bRerunMD && MASTER(cr))
1760 if (!(cr->duty & DUTY_PME))
1762 /* Tell the PME only node to finish */
1763 gmx_pme_send_finish(cr);
1768 if (ir->nstcalcenergy > 0 && !bRerunMD)
1770 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1771 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1780 pme_loadbal_done(pme_loadbal, cr, fplog, use_GPU(fr->nbv));
1783 if (shellfc && fplog)
1785 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1786 (nconverged*100.0)/step_rel);
1787 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1791 if (repl_ex_nst > 0 && MASTER(cr))
1793 print_replica_exchange_statistics(fplog, repl_ex);
1796 /* IMD cleanup, if bIMD is TRUE. */
1797 IMD_finalize(ir->bIMD, ir->imd);
1799 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);