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,2016,2017, 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.
50 #include "thread_mpi/threads.h"
52 #include "gromacs/commandline/filenm.h"
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/domdec/domdec_network.h"
55 #include "gromacs/domdec/domdec_struct.h"
56 #include "gromacs/ewald/pme.h"
57 #include "gromacs/ewald/pme-load-balancing.h"
58 #include "gromacs/fileio/trxio.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/gmxlib/nrnb.h"
61 #include "gromacs/gpu_utils/gpu_utils.h"
62 #include "gromacs/imd/imd.h"
63 #include "gromacs/listed-forces/manage-threading.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/utilities.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/math/vectypes.h"
68 #include "gromacs/mdlib/compute_io.h"
69 #include "gromacs/mdlib/constr.h"
70 #include "gromacs/mdlib/ebin.h"
71 #include "gromacs/mdlib/force.h"
72 #include "gromacs/mdlib/force_flags.h"
73 #include "gromacs/mdlib/forcerec.h"
74 #include "gromacs/mdlib/md_support.h"
75 #include "gromacs/mdlib/mdatoms.h"
76 #include "gromacs/mdlib/mdebin.h"
77 #include "gromacs/mdlib/mdoutf.h"
78 #include "gromacs/mdlib/mdrun.h"
79 #include "gromacs/mdlib/mdsetup.h"
80 #include "gromacs/mdlib/nb_verlet.h"
81 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
82 #include "gromacs/mdlib/ns.h"
83 #include "gromacs/mdlib/shellfc.h"
84 #include "gromacs/mdlib/sighandler.h"
85 #include "gromacs/mdlib/sim_util.h"
86 #include "gromacs/mdlib/simulationsignal.h"
87 #include "gromacs/mdlib/tgroup.h"
88 #include "gromacs/mdlib/trajectory_writing.h"
89 #include "gromacs/mdlib/update.h"
90 #include "gromacs/mdlib/vcm.h"
91 #include "gromacs/mdlib/vsite.h"
92 #include "gromacs/mdtypes/commrec.h"
93 #include "gromacs/mdtypes/df_history.h"
94 #include "gromacs/mdtypes/energyhistory.h"
95 #include "gromacs/mdtypes/fcdata.h"
96 #include "gromacs/mdtypes/forcerec.h"
97 #include "gromacs/mdtypes/group.h"
98 #include "gromacs/mdtypes/inputrec.h"
99 #include "gromacs/mdtypes/interaction_const.h"
100 #include "gromacs/mdtypes/md_enums.h"
101 #include "gromacs/mdtypes/mdatom.h"
102 #include "gromacs/mdtypes/state.h"
103 #include "gromacs/pbcutil/mshift.h"
104 #include "gromacs/pbcutil/pbc.h"
105 #include "gromacs/pulling/pull.h"
106 #include "gromacs/swap/swapcoords.h"
107 #include "gromacs/timing/wallcycle.h"
108 #include "gromacs/timing/walltime_accounting.h"
109 #include "gromacs/topology/atoms.h"
110 #include "gromacs/topology/idef.h"
111 #include "gromacs/topology/mtop_util.h"
112 #include "gromacs/topology/topology.h"
113 #include "gromacs/trajectory/trajectoryframe.h"
114 #include "gromacs/utility/basedefinitions.h"
115 #include "gromacs/utility/cstringutil.h"
116 #include "gromacs/utility/fatalerror.h"
117 #include "gromacs/utility/logger.h"
118 #include "gromacs/utility/real.h"
119 #include "gromacs/utility/smalloc.h"
126 #include "corewrap.h"
129 using gmx::SimulationSignaller;
131 /*! \brief Check whether bonded interactions are missing, if appropriate
133 * \param[in] fplog Log file pointer
134 * \param[in] cr Communication object
135 * \param[in] totalNumberOfBondedInteractions Result of the global reduction over the number of bonds treated in each domain
136 * \param[in] top_global Global topology for the error message
137 * \param[in] top_local Local topology for the error message
138 * \param[in] state Global state for the error message
139 * \param[inout] shouldCheckNumberOfBondedInteractions Whether we should do the check.
141 * \return Nothing, except that shouldCheckNumberOfBondedInteractions
142 * is always set to false after exit.
144 static void checkNumberOfBondedInteractions(FILE *fplog, t_commrec *cr, int totalNumberOfBondedInteractions,
145 gmx_mtop_t *top_global, gmx_localtop_t *top_local, t_state *state,
146 bool *shouldCheckNumberOfBondedInteractions)
148 if (*shouldCheckNumberOfBondedInteractions)
150 if (totalNumberOfBondedInteractions != cr->dd->nbonded_global)
152 dd_print_missing_interactions(fplog, cr, totalNumberOfBondedInteractions, top_global, top_local, state); // Does not return
154 *shouldCheckNumberOfBondedInteractions = false;
158 static void reset_all_counters(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
160 gmx_int64_t *step_rel, t_inputrec *ir,
161 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
162 gmx_walltime_accounting_t walltime_accounting,
163 struct nonbonded_verlet_t *nbv)
165 char sbuf[STEPSTRSIZE];
167 /* Reset all the counters related to performance over the run */
168 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
169 "step %s: resetting all time and cycle counters",
170 gmx_step_str(step, sbuf));
174 nbnxn_gpu_reset_timings(nbv);
178 wallcycle_stop(wcycle, ewcRUN);
179 wallcycle_reset_all(wcycle);
180 if (DOMAINDECOMP(cr))
182 reset_dd_statistics_counters(cr->dd);
185 ir->init_step += *step_rel;
186 ir->nsteps -= *step_rel;
188 wallcycle_start(wcycle, ewcRUN);
189 walltime_accounting_start(walltime_accounting);
190 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
194 \copydoc integrator_t (FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
195 int nfile, const t_filenm fnm[],
196 const gmx_output_env_t *oenv, gmx_bool bVerbose,
198 gmx_vsite_t *vsite, gmx_constr_t constr,
200 gmx::IMDOutputProvider *outputProvider,
201 t_inputrec *inputrec,
202 gmx_mtop_t *top_global, t_fcdata *fcd,
203 t_state *state_global,
205 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
208 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
209 real cpt_period, real max_hours,
212 gmx_walltime_accounting_t walltime_accounting)
214 double gmx::do_md(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
215 int nfile, const t_filenm fnm[],
216 const gmx_output_env_t *oenv, gmx_bool bVerbose,
218 gmx_vsite_t *vsite, gmx_constr_t constr,
219 int stepout, gmx::IMDOutputProvider *outputProvider,
221 gmx_mtop_t *top_global,
223 t_state *state_global,
224 energyhistory_t *energyHistory,
226 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
227 gmx_edsam_t ed, t_forcerec *fr,
228 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
229 gmx_membed_t *membed,
230 real cpt_period, real max_hours,
233 gmx_walltime_accounting_t walltime_accounting)
235 gmx_mdoutf_t outf = nullptr;
236 gmx_int64_t step, step_rel;
238 double t, t0, lam0[efptNR];
239 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
240 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD,
241 bFirstStep, startingFromCheckpoint, bInitStep, bLastStep = FALSE,
242 bBornRadii, bUsingEnsembleRestraints;
243 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
244 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
245 bForceUpdate = FALSE, bCPT;
246 gmx_bool bMasterState;
247 int force_flags, cglo_flags;
248 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
253 matrix parrinellorahmanMu, M;
255 gmx_repl_ex_t repl_ex = nullptr;
258 t_mdebin *mdebin = nullptr;
259 gmx_enerdata_t *enerd;
260 PaddedRVecVector f {};
261 gmx_global_stat_t gstat;
262 gmx_update_t *upd = nullptr;
263 t_graph *graph = nullptr;
264 gmx_groups_t *groups;
265 gmx_ekindata_t *ekind;
266 gmx_shellfc_t *shellfc;
267 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
268 gmx_bool bResetCountersHalfMaxH = FALSE;
269 gmx_bool bTemp, bPres, bTrotter;
271 rvec *cbuf = nullptr;
278 real saved_conserved_quantity = 0;
282 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
283 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
286 /* PME load balancing data for GPU kernels */
287 pme_load_balancing_t *pme_loadbal = nullptr;
288 gmx_bool bPMETune = FALSE;
289 gmx_bool bPMETunePrinting = FALSE;
292 gmx_bool bIMDstep = FALSE;
295 /* Temporary addition for FAHCORE checkpointing */
298 /* Domain decomposition could incorrectly miss a bonded
299 interaction, but checking for that requires a global
300 communication stage, which does not otherwise happen in DD
301 code. So we do that alongside the first global energy reduction
302 after a new DD is made. These variables handle whether the
303 check happens, and the result it returns. */
304 bool shouldCheckNumberOfBondedInteractions = false;
305 int totalNumberOfBondedInteractions = -1;
307 SimulationSignals signals;
308 // Most global communnication stages don't propagate mdrun
309 // signals, and will use this object to achieve that.
310 SimulationSignaller nullSignaller(nullptr, nullptr, false, false);
312 /* Check for special mdrun options */
313 bRerunMD = (Flags & MD_RERUN);
314 if (Flags & MD_RESETCOUNTERSHALFWAY)
318 /* Signal to reset the counters half the simulation steps. */
319 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
321 /* Signal to reset the counters halfway the simulation time. */
322 bResetCountersHalfMaxH = (max_hours > 0);
325 /* md-vv uses averaged full step velocities for T-control
326 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
327 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
328 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
332 /* Since we don't know if the frames read are related in any way,
333 * rebuild the neighborlist at every step.
336 ir->nstcalcenergy = 1;
340 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
341 bGStatEveryStep = (nstglobalcomm == 1);
345 ir->nstxout_compressed = 0;
347 groups = &top_global->groups;
349 if (ir->eSwapCoords != eswapNO)
351 /* Initialize ion swapping code */
352 init_swapcoords(fplog, bVerbose, ir, opt2fn_master("-swap", nfile, fnm, cr),
353 top_global, as_rvec_array(state_global->x.data()), state_global->box, &state_global->swapstate, cr, oenv, Flags);
357 init_md(fplog, cr, outputProvider, ir, oenv, &t, &t0, state_global->lambda,
358 &(state_global->fep_state), lam0,
359 nrnb, top_global, &upd,
360 nfile, fnm, &outf, &mdebin,
361 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags, wcycle);
363 clear_mat(total_vir);
365 /* Energy terms and groups */
367 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
370 /* Kinetic energy data */
372 init_ekindata(fplog, top_global, &(ir->opts), ekind);
373 /* Copy the cos acceleration to the groups struct */
374 ekind->cosacc.cos_accel = ir->cos_accel;
376 gstat = global_stat_init(ir);
378 /* Check for polarizable models and flexible constraints */
379 shellfc = init_shell_flexcon(fplog,
380 top_global, n_flexible_constraints(constr),
381 ir->nstcalcenergy, DOMAINDECOMP(cr));
383 if (shellfc && ir->nstcalcenergy != 1)
385 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);
387 if (shellfc && DOMAINDECOMP(cr))
389 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
392 if (inputrecDeform(ir))
394 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
395 set_deform_reference_box(upd,
396 deform_init_init_step_tpx,
397 deform_init_box_tpx);
398 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
402 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
403 if ((io > 2000) && MASTER(cr))
406 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
411 std::unique_ptr<t_state> stateInstance;
414 if (DOMAINDECOMP(cr))
416 top = dd_init_local_top(top_global);
418 stateInstance = std::unique_ptr<t_state>(new t_state);
419 state = stateInstance.get();
420 dd_init_local_state(cr->dd, state_global, state);
424 state_change_natoms(state_global, state_global->natoms);
425 /* We need to allocate one element extra, since we might use
426 * (unaligned) 4-wide SIMD loads to access rvec entries.
428 f.resize(state_global->natoms + 1);
429 /* Copy the pointer to the global state */
430 state = state_global;
433 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
434 &graph, mdatoms, vsite, shellfc);
436 update_realloc(upd, state->natoms);
439 /* Set up interactive MD (IMD) */
440 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, as_rvec_array(state_global->x.data()),
441 nfile, fnm, oenv, imdport, Flags);
443 if (DOMAINDECOMP(cr))
445 /* Distribute the charge groups over the nodes from the master node */
446 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
447 state_global, top_global, ir,
448 state, &f, mdatoms, top, fr,
450 nrnb, nullptr, FALSE);
451 shouldCheckNumberOfBondedInteractions = true;
452 update_realloc(upd, state->natoms);
455 update_mdatoms(mdatoms, state->lambda[efptMASS]);
457 startingFromCheckpoint = Flags & MD_STARTFROMCPT;
461 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
466 if (startingFromCheckpoint)
468 /* Update mdebin with energy history if appending to output files */
469 if (Flags & MD_APPENDFILES)
471 restore_energyhistory_from_state(mdebin, energyHistory);
475 /* We might have read an energy history from checkpoint,
476 * free the allocated memory and reset the counts.
481 /* Set the initial energy history in state by updating once */
482 update_energyhistory(energyHistory, mdebin);
485 /* Initialize constraints */
486 if (constr && !DOMAINDECOMP(cr))
488 set_constraints(constr, top, ir, mdatoms, cr);
491 if (repl_ex_nst > 0 && MASTER(cr))
493 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
494 repl_ex_nst, repl_ex_nex, repl_ex_seed);
497 /* PME tuning is only supported with PME for Coulomb. Is is not supported
498 * with only LJ PME, or for reruns.
500 bPMETune = ((Flags & MD_TUNEPME) && EEL_PME(fr->eeltype) && !bRerunMD &&
501 !(Flags & MD_REPRODUCIBLE));
504 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
505 fr->ic, fr->pmedata, use_GPU(fr->nbv),
509 if (!ir->bContinuation && !bRerunMD)
511 if (state->flags & (1 << estV))
513 /* Set the velocities of vsites, shells and frozen atoms to zero */
514 for (i = 0; i < mdatoms->homenr; i++)
516 if (mdatoms->ptype[i] == eptVSite ||
517 mdatoms->ptype[i] == eptShell)
519 clear_rvec(state->v[i]);
521 else if (mdatoms->cFREEZE)
523 for (m = 0; m < DIM; m++)
525 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
536 /* Constrain the initial coordinates and velocities */
537 do_constrain_first(fplog, constr, ir, mdatoms, state,
542 /* Construct the virtual sites for the initial configuration */
543 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, nullptr,
544 top->idef.iparams, top->idef.il,
545 fr->ePBC, fr->bMolPBC, cr, state->box);
549 if (ir->efep != efepNO)
551 /* Set free energy calculation frequency as the greatest common
552 * denominator of nstdhdl and repl_ex_nst. */
553 nstfep = ir->fepvals->nstdhdl;
556 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
560 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
564 /* Be REALLY careful about what flags you set here. You CANNOT assume
565 * this is the first step, since we might be restarting from a checkpoint,
566 * and in that case we should not do any modifications to the state.
568 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
570 if (Flags & MD_READ_EKIN)
572 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
575 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
576 | (bStopCM ? CGLO_STOPCM : 0)
577 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
578 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
579 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
581 bSumEkinhOld = FALSE;
582 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
583 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
584 constr, &nullSignaller, state->box,
585 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags
586 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
587 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
588 top_global, top, state,
589 &shouldCheckNumberOfBondedInteractions);
590 if (ir->eI == eiVVAK)
592 /* a second call to get the half step temperature initialized as well */
593 /* we do the same call as above, but turn the pressure off -- internally to
594 compute_globals, this is recognized as a velocity verlet half-step
595 kinetic energy calculation. This minimized excess variables, but
596 perhaps loses some logic?*/
598 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
599 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
600 constr, &nullSignaller, state->box,
601 nullptr, &bSumEkinhOld,
602 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
605 /* Calculate the initial half step temperature, and save the ekinh_old */
606 if (!(Flags & MD_STARTFROMCPT))
608 for (i = 0; (i < ir->opts.ngtc); i++)
610 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
615 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
616 and there is no previous step */
619 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
620 temperature control */
621 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
625 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
628 "RMS relative constraint deviation after constraining: %.2e\n",
629 constr_rmsd(constr));
631 if (EI_STATE_VELOCITY(ir->eI))
633 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
637 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
638 " input trajectory '%s'\n\n",
639 *(top_global->name), opt2fn("-rerun", nfile, fnm));
642 fprintf(stderr, "Calculated time to finish depends on nsteps from "
643 "run input file,\nwhich may not correspond to the time "
644 "needed to process input trajectory.\n\n");
650 fprintf(stderr, "starting mdrun '%s'\n",
651 *(top_global->name));
654 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
658 sprintf(tbuf, "%s", "infinite");
660 if (ir->init_step > 0)
662 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
663 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
664 gmx_step_str(ir->init_step, sbuf2),
665 ir->init_step*ir->delta_t);
669 fprintf(stderr, "%s steps, %s ps.\n",
670 gmx_step_str(ir->nsteps, sbuf), tbuf);
673 fprintf(fplog, "\n");
676 walltime_accounting_start(walltime_accounting);
677 wallcycle_start(wcycle, ewcRUN);
678 print_start(fplog, cr, walltime_accounting, "mdrun");
680 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
682 chkpt_ret = fcCheckPointParallel( cr->nodeid,
686 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
690 /***********************************************************
694 ************************************************************/
696 /* if rerunMD then read coordinates and velocities from input trajectory */
699 if (getenv("GMX_FORCE_UPDATE"))
707 bLastStep = !read_first_frame(oenv, &status,
708 opt2fn("-rerun", nfile, fnm),
709 &rerun_fr, TRX_NEED_X | TRX_READ_V);
710 if (rerun_fr.natoms != top_global->natoms)
713 "Number of atoms in trajectory (%d) does not match the "
714 "run input file (%d)\n",
715 rerun_fr.natoms, top_global->natoms);
717 if (ir->ePBC != epbcNONE)
721 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);
723 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
725 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
732 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
735 if (ir->ePBC != epbcNONE)
737 /* Set the shift vectors.
738 * Necessary here when have a static box different from the tpr box.
740 calc_shifts(rerun_fr.box, fr->shift_vec);
744 /* loop over MD steps or if rerunMD to end of input trajectory */
746 /* Skip the first Nose-Hoover integration when we get the state from tpx */
747 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
748 bSumEkinhOld = FALSE;
750 bNeedRepartition = FALSE;
751 // TODO This implementation of ensemble orientation restraints is nasty because
752 // a user can't just do multi-sim with single-sim orientation restraints.
753 bUsingEnsembleRestraints = (fcd->disres.nsystems > 1) || (cr->ms && fcd->orires.nr);
756 // Replica exchange and ensemble restraints need all
757 // simulations to remain synchronized, so they need
758 // checkpoints and stop conditions to act on the same step, so
759 // the propagation of such signals must take place between
760 // simulations, not just within simulations.
761 bool checkpointIsLocal = (repl_ex_nst <= 0) && !bUsingEnsembleRestraints;
762 bool stopConditionIsLocal = (repl_ex_nst <= 0) && !bUsingEnsembleRestraints;
763 bool resetCountersIsLocal = true;
764 signals[eglsCHKPT] = SimulationSignal(checkpointIsLocal);
765 signals[eglsSTOPCOND] = SimulationSignal(stopConditionIsLocal);
766 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
769 step = ir->init_step;
772 // TODO extract this to new multi-simulation module
773 if (MASTER(cr) && MULTISIM(cr) && (repl_ex_nst <= 0 ))
775 if (!multisim_int_all_are_equal(cr->ms, ir->nsteps))
777 GMX_LOG(mdlog.warning).appendText(
778 "Note: The number of steps is not consistent across multi simulations,\n"
779 "but we are proceeding anyway!");
781 if (!multisim_int_all_are_equal(cr->ms, ir->init_step))
783 GMX_LOG(mdlog.warning).appendText(
784 "Note: The initial step is not consistent across multi simulations,\n"
785 "but we are proceeding anyway!");
789 /* and stop now if we should */
790 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
794 /* Determine if this is a neighbor search step */
795 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
797 if (bPMETune && bNStList)
799 /* PME grid + cut-off optimization with GPUs or PME nodes */
800 pme_loadbal_do(pme_loadbal, cr,
801 (bVerbose && MASTER(cr)) ? stderr : nullptr,
809 wallcycle_start(wcycle, ewcSTEP);
815 step = rerun_fr.step;
816 step_rel = step - ir->init_step;
829 bLastStep = (step_rel == ir->nsteps);
830 t = t0 + step*ir->delta_t;
833 // TODO Refactor this, so that nstfep does not need a default value of zero
834 if (ir->efep != efepNO || ir->bSimTemp)
836 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
837 requiring different logic. */
839 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
840 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
841 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
842 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
843 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
846 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
847 do_per_step(step, repl_ex_nst));
851 update_annealing_target_temp(ir, t, upd);
856 if (!DOMAINDECOMP(cr) || MASTER(cr))
858 for (i = 0; i < state_global->natoms; i++)
860 copy_rvec(rerun_fr.x[i], state_global->x[i]);
864 for (i = 0; i < state_global->natoms; i++)
866 copy_rvec(rerun_fr.v[i], state_global->v[i]);
871 for (i = 0; i < state_global->natoms; i++)
873 clear_rvec(state_global->v[i]);
877 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
878 " Ekin, temperature and pressure are incorrect,\n"
879 " the virial will be incorrect when constraints are present.\n"
881 bRerunWarnNoV = FALSE;
885 copy_mat(rerun_fr.box, state_global->box);
886 copy_mat(state_global->box, state->box);
888 if (vsite && (Flags & MD_RERUN_VSITE))
890 if (DOMAINDECOMP(cr))
892 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
896 /* Following is necessary because the graph may get out of sync
897 * with the coordinates if we only have every N'th coordinate set
899 mk_mshift(fplog, graph, fr->ePBC, state->box, as_rvec_array(state->x.data()));
900 shift_self(graph, state->box, as_rvec_array(state->x.data()));
902 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
903 top->idef.iparams, top->idef.il,
904 fr->ePBC, fr->bMolPBC, cr, state->box);
907 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
912 /* Stop Center of Mass motion */
913 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
917 /* for rerun MD always do Neighbour Searching */
918 bNS = (bFirstStep || ir->nstlist != 0);
922 /* Determine whether or not to do Neighbour Searching */
923 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
926 /* < 0 means stop at next step, > 0 means stop at next NS step */
927 if ( (signals[eglsSTOPCOND].set < 0) ||
928 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
933 /* Determine whether or not to update the Born radii if doing GB */
934 bBornRadii = bFirstStep;
935 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
940 /* do_log triggers energy and virial calculation. Because this leads
941 * to different code paths, forces can be different. Thus for exact
942 * continuation we should avoid extra log output.
943 * Note that the || bLastStep can result in non-exact continuation
944 * beyond the last step. But we don't consider that to be an issue.
946 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
947 do_verbose = bVerbose &&
948 (step % stepout == 0 || bFirstStep || bLastStep || bRerunMD);
950 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
958 bMasterState = FALSE;
959 /* Correct the new box if it is too skewed */
960 if (inputrecDynamicBox(ir))
962 if (correct_box(fplog, step, state->box, graph))
967 if (DOMAINDECOMP(cr) && bMasterState)
969 dd_collect_state(cr->dd, state, state_global);
973 if (DOMAINDECOMP(cr))
975 /* Repartition the domain decomposition */
976 dd_partition_system(fplog, step, cr,
977 bMasterState, nstglobalcomm,
978 state_global, top_global, ir,
979 state, &f, mdatoms, top, fr,
982 do_verbose && !bPMETunePrinting);
983 shouldCheckNumberOfBondedInteractions = true;
984 update_realloc(upd, state->natoms);
988 if (MASTER(cr) && do_log)
990 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
993 if (ir->efep != efepNO)
995 update_mdatoms(mdatoms, state->lambda[efptMASS]);
998 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1001 /* We need the kinetic energy at minus the half step for determining
1002 * the full step kinetic energy and possibly for T-coupling.*/
1003 /* This may not be quite working correctly yet . . . . */
1004 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1005 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1006 constr, &nullSignaller, state->box,
1007 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1008 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1009 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1010 top_global, top, state,
1011 &shouldCheckNumberOfBondedInteractions);
1013 clear_mat(force_vir);
1015 /* We write a checkpoint at this MD step when:
1016 * either at an NS step when we signalled through gs,
1017 * or at the last step (but not when we do not want confout),
1018 * but never at the first step or with rerun.
1020 bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1021 (bLastStep && (Flags & MD_CONFOUT))) &&
1022 step > ir->init_step && !bRerunMD);
1025 signals[eglsCHKPT].set = 0;
1028 /* Determine the energy and pressure:
1029 * at nstcalcenergy steps and at energy output steps (set below).
1031 if (EI_VV(ir->eI) && (!bInitStep))
1033 /* for vv, the first half of the integration actually corresponds
1034 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1035 but the virial needs to be calculated on both the current step and the 'next' step. Future
1036 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1038 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1039 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1040 bCalcVir = bCalcEnerStep ||
1041 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1045 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1046 bCalcVir = bCalcEnerStep ||
1047 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1049 bCalcEner = bCalcEnerStep;
1051 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1053 if (do_ene || do_log || bDoReplEx)
1059 /* Do we need global communication ? */
1060 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1061 do_per_step(step, nstglobalcomm) ||
1062 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1064 force_flags = (GMX_FORCE_STATECHANGED |
1065 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1066 GMX_FORCE_ALLFORCES |
1067 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1068 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1069 (bDoFEP ? GMX_FORCE_DHDL : 0)
1074 /* Now is the time to relax the shells */
1075 relax_shell_flexcon(fplog, cr, bVerbose, step,
1076 ir, bNS, force_flags, top,
1078 state, &f, force_vir, mdatoms,
1079 nrnb, wcycle, graph, groups,
1080 shellfc, fr, bBornRadii, t, mu_tot,
1085 /* The coordinates (x) are shifted (to get whole molecules)
1087 * This is parallellized as well, and does communication too.
1088 * Check comments in sim_util.c
1090 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1091 state->box, &state->x, &state->hist,
1092 &f, force_vir, mdatoms, enerd, fcd,
1093 state->lambda, graph,
1094 fr, vsite, mu_tot, t, ed, bBornRadii,
1095 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1098 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1099 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1101 rvec *vbuf = nullptr;
1103 wallcycle_start(wcycle, ewcUPDATE);
1104 if (ir->eI == eiVV && bInitStep)
1106 /* if using velocity verlet with full time step Ekin,
1107 * take the first half step only to compute the
1108 * virial for the first step. From there,
1109 * revert back to the initial coordinates
1110 * so that the input is actually the initial step.
1112 snew(vbuf, state->natoms);
1113 copy_rvecn(as_rvec_array(state->v.data()), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1117 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1118 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1121 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1122 ekind, M, upd, etrtVELOCITY1,
1125 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1127 wallcycle_stop(wcycle, ewcUPDATE);
1128 update_constraints(fplog, step, nullptr, ir, mdatoms,
1129 state, fr->bMolPBC, graph, &f,
1130 &top->idef, shake_vir,
1131 cr, nrnb, wcycle, upd, constr,
1133 wallcycle_start(wcycle, ewcUPDATE);
1137 /* Need to unshift here if a do_force has been
1138 called in the previous step */
1139 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1141 /* if VV, compute the pressure and constraints */
1142 /* For VV2, we strictly only need this if using pressure
1143 * control, but we really would like to have accurate pressures
1145 * Think about ways around this in the future?
1146 * For now, keep this choice in comments.
1148 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1149 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1151 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1152 if (bCalcEner && ir->eI == eiVVAK)
1154 bSumEkinhOld = TRUE;
1156 /* for vv, the first half of the integration actually corresponds to the previous step.
1157 So we need information from the last step in the first half of the integration */
1158 if (bGStat || do_per_step(step-1, nstglobalcomm))
1160 wallcycle_stop(wcycle, ewcUPDATE);
1161 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1162 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1163 constr, &nullSignaller, state->box,
1164 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1165 (bGStat ? CGLO_GSTAT : 0)
1167 | (bTemp ? CGLO_TEMPERATURE : 0)
1168 | (bPres ? CGLO_PRESSURE : 0)
1169 | (bPres ? CGLO_CONSTRAINT : 0)
1170 | (bStopCM ? CGLO_STOPCM : 0)
1171 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1174 /* explanation of above:
1175 a) We compute Ekin at the full time step
1176 if 1) we are using the AveVel Ekin, and it's not the
1177 initial step, or 2) if we are using AveEkin, but need the full
1178 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1179 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1180 EkinAveVel because it's needed for the pressure */
1181 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1182 top_global, top, state,
1183 &shouldCheckNumberOfBondedInteractions);
1184 wallcycle_start(wcycle, ewcUPDATE);
1186 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1191 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1192 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1194 copy_mat(shake_vir, state->svir_prev);
1195 copy_mat(force_vir, state->fvir_prev);
1196 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1198 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1199 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1200 enerd->term[F_EKIN] = trace(ekind->ekin);
1203 else if (bExchanged)
1205 wallcycle_stop(wcycle, ewcUPDATE);
1206 /* We need the kinetic energy at minus the half step for determining
1207 * the full step kinetic energy and possibly for T-coupling.*/
1208 /* This may not be quite working correctly yet . . . . */
1209 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1210 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1211 constr, &nullSignaller, state->box,
1212 nullptr, &bSumEkinhOld,
1213 CGLO_GSTAT | CGLO_TEMPERATURE);
1214 wallcycle_start(wcycle, ewcUPDATE);
1217 /* if it's the initial step, we performed this first step just to get the constraint virial */
1218 if (ir->eI == eiVV && bInitStep)
1220 copy_rvecn(vbuf, as_rvec_array(state->v.data()), 0, state->natoms);
1223 wallcycle_stop(wcycle, ewcUPDATE);
1226 /* compute the conserved quantity */
1229 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1232 last_ekin = enerd->term[F_EKIN];
1234 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1236 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1238 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1239 if (ir->efep != efepNO && !bRerunMD)
1241 sum_dhdl(enerd, state->lambda, ir->fepvals);
1245 /* ######## END FIRST UPDATE STEP ############## */
1246 /* ######## If doing VV, we now have v(dt) ###### */
1249 /* perform extended ensemble sampling in lambda - we don't
1250 actually move to the new state before outputting
1251 statistics, but if performing simulated tempering, we
1252 do update the velocities and the tau_t. */
1254 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, as_rvec_array(state->v.data()), mdatoms);
1255 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1256 copy_df_history(state_global->dfhist, state->dfhist);
1259 /* Now we have the energies and forces corresponding to the
1260 * coordinates at time t. We must output all of this before
1263 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1264 ir, state, state_global, energyHistory,
1266 outf, mdebin, ekind, &f,
1268 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1270 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1271 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, as_rvec_array(state->x.data()), ir, t, wcycle);
1273 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1274 if (startingFromCheckpoint && bTrotter)
1276 copy_mat(state->svir_prev, shake_vir);
1277 copy_mat(state->fvir_prev, force_vir);
1280 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1282 /* Check whether everything is still allright */
1283 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1289 int nsteps_stop = -1;
1291 /* this just makes signals[].sig compatible with the hack
1292 of sending signals around by MPI_Reduce together with
1294 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1296 signals[eglsSTOPCOND].sig = 1;
1297 nsteps_stop = std::max(ir->nstlist, 2*nstglobalcomm);
1299 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1301 signals[eglsSTOPCOND].sig = -1;
1302 nsteps_stop = nstglobalcomm + 1;
1307 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1308 gmx_get_signal_name(), nsteps_stop);
1312 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1313 gmx_get_signal_name(), nsteps_stop);
1315 handled_stop_condition = (int)gmx_get_stop_condition();
1317 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1318 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1319 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1321 /* Signal to terminate the run */
1322 signals[eglsSTOPCOND].sig = 1;
1325 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1327 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1330 if (bResetCountersHalfMaxH && MASTER(cr) &&
1331 elapsed_time > max_hours*60.0*60.0*0.495)
1333 /* Set flag that will communicate the signal to all ranks in the simulation */
1334 signals[eglsRESETCOUNTERS].sig = 1;
1337 /* In parallel we only have to check for checkpointing in steps
1338 * where we do global communication,
1339 * otherwise the other nodes don't know.
1341 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1344 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1345 signals[eglsCHKPT].set == 0)
1347 signals[eglsCHKPT].sig = 1;
1350 /* ######### START SECOND UPDATE STEP ################# */
1352 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1355 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1357 gmx_bool bIfRandomize;
1358 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1359 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1360 if (constr && bIfRandomize)
1362 update_constraints(fplog, step, nullptr, ir, mdatoms,
1363 state, fr->bMolPBC, graph, &f,
1364 &top->idef, tmp_vir,
1365 cr, nrnb, wcycle, upd, constr,
1369 /* Box is changed in update() when we do pressure coupling,
1370 * but we should still use the old box for energy corrections and when
1371 * writing it to the energy file, so it matches the trajectory files for
1372 * the same timestep above. Make a copy in a separate array.
1374 copy_mat(state->box, lastbox);
1378 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1380 wallcycle_start(wcycle, ewcUPDATE);
1381 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1384 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1385 /* We can only do Berendsen coupling after we have summed
1386 * the kinetic energy or virial. Since the happens
1387 * in global_state after update, we should only do it at
1388 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1393 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1394 update_pcouple_before_coordinates(fplog, step, ir, state,
1395 parrinellorahmanMu, M,
1401 /* velocity half-step update */
1402 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1403 ekind, M, upd, etrtVELOCITY2,
1407 /* Above, initialize just copies ekinh into ekin,
1408 * it doesn't copy position (for VV),
1409 * and entire integrator for MD.
1412 if (ir->eI == eiVVAK)
1414 /* We probably only need md->homenr, not state->natoms */
1415 if (state->natoms > cbuf_nalloc)
1417 cbuf_nalloc = state->natoms;
1418 srenew(cbuf, cbuf_nalloc);
1420 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1423 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1424 ekind, M, upd, etrtPOSITION, cr, constr);
1425 wallcycle_stop(wcycle, ewcUPDATE);
1427 update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1428 fr->bMolPBC, graph, &f,
1429 &top->idef, shake_vir,
1430 cr, nrnb, wcycle, upd, constr,
1433 if (ir->eI == eiVVAK)
1435 /* erase F_EKIN and F_TEMP here? */
1436 /* just compute the kinetic energy at the half step to perform a trotter step */
1437 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1438 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1439 constr, &nullSignaller, lastbox,
1440 nullptr, &bSumEkinhOld,
1441 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1443 wallcycle_start(wcycle, ewcUPDATE);
1444 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1445 /* now we know the scaling, we can compute the positions again again */
1446 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1448 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1449 ekind, M, upd, etrtPOSITION, cr, constr);
1450 wallcycle_stop(wcycle, ewcUPDATE);
1452 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1453 /* are the small terms in the shake_vir here due
1454 * to numerical errors, or are they important
1455 * physically? I'm thinking they are just errors, but not completely sure.
1456 * For now, will call without actually constraining, constr=NULL*/
1457 update_constraints(fplog, step, nullptr, ir, mdatoms,
1458 state, fr->bMolPBC, graph, &f,
1459 &top->idef, tmp_vir,
1460 cr, nrnb, wcycle, upd, nullptr,
1465 /* this factor or 2 correction is necessary
1466 because half of the constraint force is removed
1467 in the vv step, so we have to double it. See
1468 the Redmine issue #1255. It is not yet clear
1469 if the factor of 2 is exact, or just a very
1470 good approximation, and this will be
1471 investigated. The next step is to see if this
1472 can be done adding a dhdl contribution from the
1473 rattle step, but this is somewhat more
1474 complicated with the current code. Will be
1475 investigated, hopefully for 4.6.3. However,
1476 this current solution is much better than
1477 having it completely wrong.
1479 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1483 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1488 /* Need to unshift here */
1489 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1492 if (vsite != nullptr)
1494 wallcycle_start(wcycle, ewcVSITECONSTR);
1495 if (graph != nullptr)
1497 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1499 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
1500 top->idef.iparams, top->idef.il,
1501 fr->ePBC, fr->bMolPBC, cr, state->box);
1503 if (graph != nullptr)
1505 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1507 wallcycle_stop(wcycle, ewcVSITECONSTR);
1510 /* ############## IF NOT VV, Calculate globals HERE ############ */
1511 /* With Leap-Frog we can skip compute_globals at
1512 * non-communication steps, but we need to calculate
1513 * the kinetic energy one step before communication.
1516 // Organize to do inter-simulation signalling on steps if
1517 // and when algorithms require it.
1518 bool doInterSimSignal = (!bFirstStep && bDoReplEx) || bUsingEnsembleRestraints;
1520 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1522 // Since we're already communicating at this step, we
1523 // can propagate intra-simulation signals. Note that
1524 // check_nstglobalcomm has the responsibility for
1525 // choosing the value of nstglobalcomm that is one way
1526 // bGStat becomes true, so we can't get into a
1527 // situation where e.g. checkpointing can't be
1529 bool doIntraSimSignal = true;
1530 SimulationSignaller signaller(&signals, cr, doInterSimSignal, doIntraSimSignal);
1532 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1533 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1536 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1537 (bGStat ? CGLO_GSTAT : 0)
1538 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1539 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1540 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1541 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1543 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1545 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1546 top_global, top, state,
1547 &shouldCheckNumberOfBondedInteractions);
1551 /* ############# END CALC EKIN AND PRESSURE ################# */
1553 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1554 the virial that should probably be addressed eventually. state->veta has better properies,
1555 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1556 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1558 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1560 /* Sum up the foreign energy and dhdl terms for md and sd.
1561 Currently done every step so that dhdl is correct in the .edr */
1562 sum_dhdl(enerd, state->lambda, ir->fepvals);
1565 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1566 pres, parrinellorahmanMu,
1569 /* ################# END UPDATE STEP 2 ################# */
1570 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1572 /* The coordinates (x) were unshifted in update */
1575 /* We will not sum ekinh_old,
1576 * so signal that we still have to do it.
1578 bSumEkinhOld = TRUE;
1581 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1583 /* use the directly determined last velocity, not actually the averaged half steps */
1584 if (bTrotter && ir->eI == eiVV)
1586 enerd->term[F_EKIN] = last_ekin;
1588 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1592 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1596 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1598 /* ######### END PREPARING EDR OUTPUT ########### */
1603 if (fplog && do_log && bDoExpanded)
1605 /* only needed if doing expanded ensemble */
1606 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1607 state_global->dfhist, state->fep_state, ir->nstlog, step);
1611 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1612 t, mdatoms->tmass, enerd, state,
1613 ir->fepvals, ir->expandedvals, lastbox,
1614 shake_vir, force_vir, total_vir, pres,
1615 ekind, mu_tot, constr);
1619 upd_mdebin_step(mdebin);
1622 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1623 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1625 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1627 eprNORMAL, mdebin, fcd, groups, &(ir->opts));
1631 pull_print_output(ir->pull_work, step, t);
1634 if (do_per_step(step, ir->nstlog))
1636 if (fflush(fplog) != 0)
1638 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1644 /* Have to do this part _after_ outputting the logfile and the edr file */
1645 /* Gets written into the state at the beginning of next loop*/
1646 state->fep_state = lamnew;
1648 /* Print the remaining wall clock time for the run */
1649 if (MULTIMASTER(cr) &&
1650 (do_verbose || gmx_got_usr_signal()) &&
1655 fprintf(stderr, "\n");
1657 print_time(stderr, walltime_accounting, step, ir, cr);
1660 /* Ion/water position swapping.
1661 * Not done in last step since trajectory writing happens before this call
1662 * in the MD loop and exchanges would be lost anyway. */
1663 bNeedRepartition = FALSE;
1664 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1665 do_per_step(step, ir->swap->nstswap))
1667 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1668 bRerunMD ? rerun_fr.x : as_rvec_array(state->x.data()),
1669 bRerunMD ? rerun_fr.box : state->box,
1670 MASTER(cr) && bVerbose, bRerunMD);
1672 if (bNeedRepartition && DOMAINDECOMP(cr))
1674 dd_collect_state(cr->dd, state, state_global);
1678 /* Replica exchange */
1682 bExchanged = replica_exchange(fplog, cr, repl_ex,
1683 state_global, enerd,
1687 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1689 dd_partition_system(fplog, step, cr, TRUE, 1,
1690 state_global, top_global, ir,
1691 state, &f, mdatoms, top, fr,
1693 nrnb, wcycle, FALSE);
1694 shouldCheckNumberOfBondedInteractions = true;
1695 update_realloc(upd, state->natoms);
1700 startingFromCheckpoint = FALSE;
1702 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1703 /* With all integrators, except VV, we need to retain the pressure
1704 * at the current step for coupling at the next step.
1706 if ((state->flags & (1<<estPRES_PREV)) &&
1708 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1710 /* Store the pressure in t_state for pressure coupling
1711 * at the next MD step.
1713 copy_mat(pres, state->pres_prev);
1716 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1718 if ( (membed != nullptr) && (!bLastStep) )
1720 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1727 /* read next frame from input trajectory */
1728 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1733 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1737 cycles = wallcycle_stop(wcycle, ewcSTEP);
1738 if (DOMAINDECOMP(cr) && wcycle)
1740 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1743 if (!bRerunMD || !rerun_fr.bStep)
1745 /* increase the MD step number */
1750 /* TODO make a counter-reset module */
1751 /* If it is time to reset counters, set a flag that remains
1752 true until counters actually get reset */
1753 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1754 signals[eglsRESETCOUNTERS].set != 0)
1756 if (pme_loadbal_is_active(pme_loadbal))
1758 /* Do not permit counter reset while PME load
1759 * balancing is active. The only purpose for resetting
1760 * counters is to measure reliable performance data,
1761 * and that can't be done before balancing
1764 * TODO consider fixing this by delaying the reset
1765 * until after load balancing completes,
1766 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1767 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1768 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1769 "resetting counters later in the run, e.g. with gmx "
1770 "mdrun -resetstep.", step);
1772 reset_all_counters(fplog, mdlog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1773 use_GPU(fr->nbv) ? fr->nbv : nullptr);
1774 wcycle_set_reset_counters(wcycle, -1);
1775 if (!(cr->duty & DUTY_PME))
1777 /* Tell our PME node to reset its counters */
1778 gmx_pme_send_resetcounters(cr, step);
1780 /* Correct max_hours for the elapsed time */
1781 max_hours -= elapsed_time/(60.0*60.0);
1782 /* If mdrun -maxh -resethway was active, it can only trigger once */
1783 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1784 /* Reset can only happen once, so clear the triggering flag. */
1785 signals[eglsRESETCOUNTERS].set = 0;
1788 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1789 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1792 /* End of main MD loop */
1794 /* Closing TNG files can include compressing data. Therefore it is good to do that
1795 * before stopping the time measurements. */
1796 mdoutf_tng_close(outf);
1798 /* Stop measuring walltime */
1799 walltime_accounting_end(walltime_accounting);
1801 if (bRerunMD && MASTER(cr))
1806 if (!(cr->duty & DUTY_PME))
1808 /* Tell the PME only node to finish */
1809 gmx_pme_send_finish(cr);
1814 if (ir->nstcalcenergy > 0 && !bRerunMD)
1816 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1817 eprAVER, mdebin, fcd, groups, &(ir->opts));
1825 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
1828 done_shellfc(fplog, shellfc, step_rel);
1830 if (repl_ex_nst > 0 && MASTER(cr))
1832 print_replica_exchange_statistics(fplog, repl_ex);
1835 // Clean up swapcoords
1836 if (ir->eSwapCoords != eswapNO)
1838 finish_swapcoords(ir->swap);
1841 /* IMD cleanup, if bIMD is TRUE. */
1842 IMD_finalize(ir->bIMD, ir->imd);
1844 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1845 if (step_rel >= wcycle_get_reset_counters(wcycle) &&
1846 signals[eglsRESETCOUNTERS].set == 0 &&
1847 !bResetCountersHalfMaxH)
1849 walltime_accounting_set_valid_finish(walltime_accounting);