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,2018, 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.
51 #include "thread_mpi/threads.h"
53 #include "gromacs/awh/awh.h"
54 #include "gromacs/commandline/filenm.h"
55 #include "gromacs/compat/make_unique.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/essentialdynamics/edsam.h"
60 #include "gromacs/ewald/pme.h"
61 #include "gromacs/ewald/pme-load-balancing.h"
62 #include "gromacs/fileio/trxio.h"
63 #include "gromacs/gmxlib/network.h"
64 #include "gromacs/gmxlib/nrnb.h"
65 #include "gromacs/gpu_utils/gpu_utils.h"
66 #include "gromacs/imd/imd.h"
67 #include "gromacs/listed-forces/manage-threading.h"
68 #include "gromacs/math/functions.h"
69 #include "gromacs/math/utilities.h"
70 #include "gromacs/math/vec.h"
71 #include "gromacs/math/vectypes.h"
72 #include "gromacs/mdlib/compute_io.h"
73 #include "gromacs/mdlib/constr.h"
74 #include "gromacs/mdlib/ebin.h"
75 #include "gromacs/mdlib/expanded.h"
76 #include "gromacs/mdlib/force.h"
77 #include "gromacs/mdlib/force_flags.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/md_support.h"
80 #include "gromacs/mdlib/mdatoms.h"
81 #include "gromacs/mdlib/mdebin.h"
82 #include "gromacs/mdlib/mdoutf.h"
83 #include "gromacs/mdlib/mdrun.h"
84 #include "gromacs/mdlib/mdsetup.h"
85 #include "gromacs/mdlib/nb_verlet.h"
86 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
87 #include "gromacs/mdlib/ns.h"
88 #include "gromacs/mdlib/shellfc.h"
89 #include "gromacs/mdlib/sighandler.h"
90 #include "gromacs/mdlib/sim_util.h"
91 #include "gromacs/mdlib/simulationsignal.h"
92 #include "gromacs/mdlib/tgroup.h"
93 #include "gromacs/mdlib/trajectory_writing.h"
94 #include "gromacs/mdlib/update.h"
95 #include "gromacs/mdlib/vcm.h"
96 #include "gromacs/mdlib/vsite.h"
97 #include "gromacs/mdtypes/awh-history.h"
98 #include "gromacs/mdtypes/commrec.h"
99 #include "gromacs/mdtypes/df_history.h"
100 #include "gromacs/mdtypes/energyhistory.h"
101 #include "gromacs/mdtypes/fcdata.h"
102 #include "gromacs/mdtypes/forcerec.h"
103 #include "gromacs/mdtypes/group.h"
104 #include "gromacs/mdtypes/inputrec.h"
105 #include "gromacs/mdtypes/interaction_const.h"
106 #include "gromacs/mdtypes/md_enums.h"
107 #include "gromacs/mdtypes/mdatom.h"
108 #include "gromacs/mdtypes/observableshistory.h"
109 #include "gromacs/mdtypes/state.h"
110 #include "gromacs/pbcutil/mshift.h"
111 #include "gromacs/pbcutil/pbc.h"
112 #include "gromacs/pulling/pull.h"
113 #include "gromacs/swap/swapcoords.h"
114 #include "gromacs/timing/wallcycle.h"
115 #include "gromacs/timing/walltime_accounting.h"
116 #include "gromacs/topology/atoms.h"
117 #include "gromacs/topology/idef.h"
118 #include "gromacs/topology/mtop_util.h"
119 #include "gromacs/topology/topology.h"
120 #include "gromacs/trajectory/trajectoryframe.h"
121 #include "gromacs/utility/basedefinitions.h"
122 #include "gromacs/utility/cstringutil.h"
123 #include "gromacs/utility/fatalerror.h"
124 #include "gromacs/utility/logger.h"
125 #include "gromacs/utility/real.h"
126 #include "gromacs/utility/smalloc.h"
133 #include "corewrap.h"
136 using gmx::SimulationSignaller;
138 /*! \brief Check whether bonded interactions are missing, if appropriate
140 * \param[in] fplog Log file pointer
141 * \param[in] cr Communication object
142 * \param[in] totalNumberOfBondedInteractions Result of the global reduction over the number of bonds treated in each domain
143 * \param[in] top_global Global topology for the error message
144 * \param[in] top_local Local topology for the error message
145 * \param[in] state Global state for the error message
146 * \param[inout] shouldCheckNumberOfBondedInteractions Whether we should do the check.
148 * \return Nothing, except that shouldCheckNumberOfBondedInteractions
149 * is always set to false after exit.
151 static void checkNumberOfBondedInteractions(FILE *fplog, t_commrec *cr, int totalNumberOfBondedInteractions,
152 gmx_mtop_t *top_global, gmx_localtop_t *top_local, t_state *state,
153 bool *shouldCheckNumberOfBondedInteractions)
155 if (*shouldCheckNumberOfBondedInteractions)
157 if (totalNumberOfBondedInteractions != cr->dd->nbonded_global)
159 dd_print_missing_interactions(fplog, cr, totalNumberOfBondedInteractions, top_global, top_local, state); // Does not return
161 *shouldCheckNumberOfBondedInteractions = false;
165 static void reset_all_counters(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
167 gmx_int64_t *step_rel, t_inputrec *ir,
168 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
169 gmx_walltime_accounting_t walltime_accounting,
170 struct nonbonded_verlet_t *nbv,
171 struct gmx_pme_t *pme)
173 char sbuf[STEPSTRSIZE];
175 /* Reset all the counters related to performance over the run */
176 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
177 "step %s: resetting all time and cycle counters",
178 gmx_step_str(step, sbuf));
182 nbnxn_gpu_reset_timings(nbv);
185 if (pme_gpu_task_enabled(pme))
187 pme_gpu_reset_timings(pme);
190 if (use_GPU(nbv) || pme_gpu_task_enabled(pme))
195 wallcycle_stop(wcycle, ewcRUN);
196 wallcycle_reset_all(wcycle);
197 if (DOMAINDECOMP(cr))
199 reset_dd_statistics_counters(cr->dd);
202 ir->init_step += *step_rel;
203 ir->nsteps -= *step_rel;
205 wallcycle_start(wcycle, ewcRUN);
206 walltime_accounting_start(walltime_accounting);
207 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
210 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
212 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
213 * \param[in,out] globalState The global state container
214 * \param[in] constructVsites When true, vsite coordinates are constructed
215 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
216 * \param[in] idef Topology parameters, used for constructing vsites
217 * \param[in] timeStep Time step, used for constructing vsites
218 * \param[in] forceRec Force record, used for constructing vsites
219 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
220 * \param[in,out] warnWhenNoV When true, issue a warning when no velocities are present in \p rerunFrame; is set to false when a warning was issued
222 static void prepareRerunState(const t_trxframe &rerunFrame,
223 t_state *globalState,
224 bool constructVsites,
225 const gmx_vsite_t *vsite,
228 const t_forcerec &forceRec,
230 gmx_bool *warnWhenNoV)
232 for (int i = 0; i < globalState->natoms; i++)
234 copy_rvec(rerunFrame.x[i], globalState->x[i]);
238 for (int i = 0; i < globalState->natoms; i++)
240 copy_rvec(rerunFrame.v[i], globalState->v[i]);
245 for (int i = 0; i < globalState->natoms; i++)
247 clear_rvec(globalState->v[i]);
251 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
252 " Ekin, temperature and pressure are incorrect,\n"
253 " the virial will be incorrect when constraints are present.\n"
255 *warnWhenNoV = FALSE;
258 copy_mat(rerunFrame.box, globalState->box);
262 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
266 /* Following is necessary because the graph may get out of sync
267 * with the coordinates if we only have every N'th coordinate set
269 mk_mshift(nullptr, graph, forceRec.ePBC, globalState->box, as_rvec_array(globalState->x.data()));
270 shift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
272 construct_vsites(vsite, as_rvec_array(globalState->x.data()), timeStep, as_rvec_array(globalState->v.data()),
273 idef.iparams, idef.il,
274 forceRec.ePBC, forceRec.bMolPBC, nullptr, globalState->box);
277 unshift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
283 \copydoc integrator_t (FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
284 int nfile, const t_filenm fnm[],
285 const gmx_output_env_t *oenv,
286 const MdrunOptions &mdrunOptions,
287 gmx_vsite_t *vsite, gmx_constr_t constr,
288 gmx::IMDOutputProvider *outputProvider,
289 t_inputrec *inputrec,
290 gmx_mtop_t *top_global, t_fcdata *fcd,
291 t_state *state_global,
293 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
295 const ReplicaExchangeParameters &replExParams,
296 gmx_membed_t *membed,
297 gmx_walltime_accounting_t walltime_accounting)
299 double gmx::do_md(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
300 int nfile, const t_filenm fnm[],
301 const gmx_output_env_t *oenv,
302 const MdrunOptions &mdrunOptions,
303 gmx_vsite_t *vsite, gmx_constr_t constr,
304 gmx::IMDOutputProvider *outputProvider,
306 gmx_mtop_t *top_global,
308 t_state *state_global,
309 ObservablesHistory *observablesHistory,
310 gmx::MDAtoms *mdAtoms,
311 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
313 const ReplicaExchangeParameters &replExParams,
314 gmx_membed_t *membed,
315 gmx_walltime_accounting_t walltime_accounting)
317 gmx_mdoutf_t outf = nullptr;
318 gmx_int64_t step, step_rel;
320 double t, t0, lam0[efptNR];
321 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
322 gmx_bool bNS, bNStList, bSimAnn, bStopCM,
323 bFirstStep, bInitStep, bLastStep = FALSE,
324 bBornRadii, bUsingEnsembleRestraints;
325 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
326 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
327 bForceUpdate = FALSE, bCPT;
328 gmx_bool bMasterState;
329 int force_flags, cglo_flags;
330 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
335 matrix parrinellorahmanMu, M;
337 gmx_repl_ex_t repl_ex = nullptr;
340 t_mdebin *mdebin = nullptr;
341 gmx_enerdata_t *enerd;
342 PaddedRVecVector f {};
343 gmx_global_stat_t gstat;
344 gmx_update_t *upd = nullptr;
345 t_graph *graph = nullptr;
346 gmx_groups_t *groups;
347 gmx_ekindata_t *ekind;
348 gmx_shellfc_t *shellfc;
349 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
350 gmx_bool bResetCountersHalfMaxH = FALSE;
351 gmx_bool bTemp, bPres, bTrotter;
353 rvec *cbuf = nullptr;
360 real saved_conserved_quantity = 0;
364 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
365 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
368 /* PME load balancing data for GPU kernels */
369 pme_load_balancing_t *pme_loadbal = nullptr;
370 gmx_bool bPMETune = FALSE;
371 gmx_bool bPMETunePrinting = FALSE;
374 gmx_bool bIMDstep = FALSE;
377 /* Temporary addition for FAHCORE checkpointing */
380 /* Domain decomposition could incorrectly miss a bonded
381 interaction, but checking for that requires a global
382 communication stage, which does not otherwise happen in DD
383 code. So we do that alongside the first global energy reduction
384 after a new DD is made. These variables handle whether the
385 check happens, and the result it returns. */
386 bool shouldCheckNumberOfBondedInteractions = false;
387 int totalNumberOfBondedInteractions = -1;
389 SimulationSignals signals;
390 // Most global communnication stages don't propagate mdrun
391 // signals, and will use this object to achieve that.
392 SimulationSignaller nullSignaller(nullptr, nullptr, false, false);
394 if (!mdrunOptions.writeConfout)
396 // This is on by default, and the main known use case for
397 // turning it off is for convenience in benchmarking, which is
398 // something that should not show up in the general user
400 GMX_LOG(mdlog.info).asParagraph().
401 appendText("The -noconfout functionality is deprecated, and may be removed in a future version.");
403 if (mdrunOptions.timingOptions.resetHalfway)
405 GMX_LOG(mdlog.info).asParagraph().
406 appendText("The -resethway functionality is deprecated, and may be removed in a future version.");
409 /* Signal to reset the counters half the simulation steps. */
410 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
412 /* Signal to reset the counters halfway the simulation time. */
413 bResetCountersHalfMaxH = (mdrunOptions.maximumHoursToRun > 0);
416 /* md-vv uses averaged full step velocities for T-control
417 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
418 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
419 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
421 const gmx_bool bRerunMD = mdrunOptions.rerun;
422 int nstglobalcomm = mdrunOptions.globalCommunicationInterval;
426 /* Since we don't know if the frames read are related in any way,
427 * rebuild the neighborlist at every step.
430 ir->nstcalcenergy = 1;
434 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
435 bGStatEveryStep = (nstglobalcomm == 1);
439 ir->nstxout_compressed = 0;
441 groups = &top_global->groups;
443 gmx_edsam *ed = nullptr;
444 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
446 /* Initialize essential dynamics sampling */
447 ed = init_edsam(opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm),
450 state_global, observablesHistory,
451 oenv, mdrunOptions.continuationOptions.appendFiles);
454 if (ir->eSwapCoords != eswapNO)
456 /* Initialize ion swapping code */
457 init_swapcoords(fplog, ir, opt2fn_master("-swap", nfile, fnm, cr),
459 state_global, observablesHistory,
460 cr, oenv, mdrunOptions);
464 init_md(fplog, cr, outputProvider, ir, oenv, mdrunOptions,
465 &t, &t0, state_global, lam0,
466 nrnb, top_global, &upd,
467 nfile, fnm, &outf, &mdebin,
468 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, wcycle);
470 clear_mat(total_vir);
472 /* Energy terms and groups */
474 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
477 /* Kinetic energy data */
479 init_ekindata(fplog, top_global, &(ir->opts), ekind);
480 /* Copy the cos acceleration to the groups struct */
481 ekind->cosacc.cos_accel = ir->cos_accel;
483 gstat = global_stat_init(ir);
485 /* Check for polarizable models and flexible constraints */
486 shellfc = init_shell_flexcon(fplog,
487 top_global, n_flexible_constraints(constr),
488 ir->nstcalcenergy, DOMAINDECOMP(cr));
490 if (shellfc && ir->nstcalcenergy != 1)
492 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);
494 if (shellfc && DOMAINDECOMP(cr))
496 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
498 if (shellfc && ir->bDoAwh)
500 gmx_fatal(FARGS, "AWH biasing does not support shell particles.");
503 if (inputrecDeform(ir))
505 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
506 set_deform_reference_box(upd,
507 deform_init_init_step_tpx,
508 deform_init_box_tpx);
509 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
513 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
514 if ((io > 2000) && MASTER(cr))
517 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
522 // Local state only becomes valid now.
523 std::unique_ptr<t_state> stateInstance;
526 if (DOMAINDECOMP(cr))
528 top = dd_init_local_top(top_global);
530 stateInstance = gmx::compat::make_unique<t_state>();
531 state = stateInstance.get();
532 dd_init_local_state(cr->dd, state_global, state);
536 state_change_natoms(state_global, state_global->natoms);
537 /* We need to allocate one element extra, since we might use
538 * (unaligned) 4-wide SIMD loads to access rvec entries.
540 f.resize(gmx::paddedRVecVectorSize(state_global->natoms));
541 /* Copy the pointer to the global state */
542 state = state_global;
545 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
546 &graph, mdAtoms, vsite, shellfc);
548 update_realloc(upd, state->natoms);
551 /* Set up interactive MD (IMD) */
552 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
553 nfile, fnm, oenv, mdrunOptions);
555 if (DOMAINDECOMP(cr))
557 /* Distribute the charge groups over the nodes from the master node */
558 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
559 state_global, top_global, ir,
560 state, &f, mdAtoms, top, fr,
562 nrnb, nullptr, FALSE);
563 shouldCheckNumberOfBondedInteractions = true;
564 update_realloc(upd, state->natoms);
567 auto mdatoms = mdAtoms->mdatoms();
569 // NOTE: The global state is no longer used at this point.
570 // But state_global is still used as temporary storage space for writing
571 // the global state to file and potentially for replica exchange.
572 // (Global topology should persist.)
574 update_mdatoms(mdatoms, state->lambda[efptMASS]);
576 const ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
577 bool startingFromCheckpoint = continuationOptions.startedFromCheckpoint;
581 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
586 if (startingFromCheckpoint)
588 /* Update mdebin with energy history if appending to output files */
589 if (continuationOptions.appendFiles)
591 restore_energyhistory_from_state(mdebin, observablesHistory->energyHistory.get());
593 else if (observablesHistory->energyHistory.get() != nullptr)
595 /* We might have read an energy history from checkpoint.
596 * As we are not appending, we want to restart the statistics.
597 * Free the allocated memory and reset the counts.
599 observablesHistory->energyHistory = {};
602 if (observablesHistory->energyHistory.get() == nullptr)
604 observablesHistory->energyHistory = std::unique_ptr<energyhistory_t>(new energyhistory_t {});
606 /* Set the initial energy history in state by updating once */
607 update_energyhistory(observablesHistory->energyHistory.get(), mdebin);
610 /* Initialize constraints */
611 if (constr && !DOMAINDECOMP(cr))
613 set_constraints(constr, top, ir, mdatoms, cr);
616 /* Initialize AWH and restore state from history in checkpoint if needed. */
619 ir->awh = new gmx::Awh(fplog, *ir, cr, *ir->awhParams, opt2fn("-awh", nfile, fnm), ir->pull_work);
621 if (startingFromCheckpoint)
623 /* Restore the AWH history read from checkpoint */
624 ir->awh->restoreStateFromHistory(MASTER(cr) ? state_global->awhHistory.get() : nullptr);
628 /* Initialize the AWH history here */
629 state_global->awhHistory = ir->awh->initHistoryFromState();
633 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
634 if (useReplicaExchange && MASTER(cr))
636 repl_ex = init_replica_exchange(fplog, cr->ms, top_global->natoms, ir,
640 /* PME tuning is only supported in the Verlet scheme, with PME for
641 * Coulomb. It is not supported with only LJ PME, or for
643 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !bRerunMD &&
644 !mdrunOptions.reproducible && ir->cutoff_scheme != ecutsGROUP);
647 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
648 fr->ic, fr->nbv->listParams.get(), fr->pmedata, use_GPU(fr->nbv),
652 if (!ir->bContinuation && !bRerunMD)
654 if (state->flags & (1 << estV))
656 /* Set the velocities of vsites, shells and frozen atoms to zero */
657 for (i = 0; i < mdatoms->homenr; i++)
659 if (mdatoms->ptype[i] == eptVSite ||
660 mdatoms->ptype[i] == eptShell)
662 clear_rvec(state->v[i]);
664 else if (mdatoms->cFREEZE)
666 for (m = 0; m < DIM; m++)
668 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
679 /* Constrain the initial coordinates and velocities */
680 do_constrain_first(fplog, constr, ir, mdatoms, state,
685 /* Construct the virtual sites for the initial configuration */
686 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, nullptr,
687 top->idef.iparams, top->idef.il,
688 fr->ePBC, fr->bMolPBC, cr, state->box);
692 if (ir->efep != efepNO)
694 /* Set free energy calculation frequency as the greatest common
695 * denominator of nstdhdl and repl_ex_nst. */
696 nstfep = ir->fepvals->nstdhdl;
699 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
701 if (useReplicaExchange)
703 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
707 /* Be REALLY careful about what flags you set here. You CANNOT assume
708 * this is the first step, since we might be restarting from a checkpoint,
709 * and in that case we should not do any modifications to the state.
711 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
713 if (continuationOptions.haveReadEkin)
715 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
718 cglo_flags = (CGLO_INITIALIZATION | CGLO_TEMPERATURE | CGLO_GSTAT
719 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
720 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
721 | (continuationOptions.haveReadEkin ? CGLO_READEKIN : 0));
723 bSumEkinhOld = FALSE;
724 /* To minimize communication, compute_globals computes the COM velocity
725 * and the kinetic energy for the velocities without COM motion removed.
726 * Thus to get the kinetic energy without the COM contribution, we need
727 * to call compute_globals twice.
729 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
731 int cglo_flags_iteration = cglo_flags;
732 if (bStopCM && cgloIteration == 0)
734 cglo_flags_iteration |= CGLO_STOPCM;
735 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
737 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
738 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
739 constr, &nullSignaller, state->box,
740 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags_iteration
741 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
743 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
744 top_global, top, state,
745 &shouldCheckNumberOfBondedInteractions);
746 if (ir->eI == eiVVAK)
748 /* a second call to get the half step temperature initialized as well */
749 /* we do the same call as above, but turn the pressure off -- internally to
750 compute_globals, this is recognized as a velocity verlet half-step
751 kinetic energy calculation. This minimized excess variables, but
752 perhaps loses some logic?*/
754 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
755 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
756 constr, &nullSignaller, state->box,
757 nullptr, &bSumEkinhOld,
758 cglo_flags & ~CGLO_PRESSURE);
761 /* Calculate the initial half step temperature, and save the ekinh_old */
762 if (!continuationOptions.startedFromCheckpoint)
764 for (i = 0; (i < ir->opts.ngtc); i++)
766 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
770 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
771 temperature control */
772 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
776 if (!ir->bContinuation)
778 if (constr && ir->eConstrAlg == econtLINCS)
781 "RMS relative constraint deviation after constraining: %.2e\n",
782 constr_rmsd(constr));
784 if (EI_STATE_VELOCITY(ir->eI))
786 real temp = enerd->term[F_TEMP];
789 /* Result of Ekin averaged over velocities of -half
790 * and +half step, while we only have -half step here.
794 fprintf(fplog, "Initial temperature: %g K\n", temp);
800 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
801 " input trajectory '%s'\n\n",
802 *(top_global->name), opt2fn("-rerun", nfile, fnm));
803 if (mdrunOptions.verbose)
805 fprintf(stderr, "Calculated time to finish depends on nsteps from "
806 "run input file,\nwhich may not correspond to the time "
807 "needed to process input trajectory.\n\n");
813 fprintf(stderr, "starting mdrun '%s'\n",
814 *(top_global->name));
817 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
821 sprintf(tbuf, "%s", "infinite");
823 if (ir->init_step > 0)
825 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
826 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
827 gmx_step_str(ir->init_step, sbuf2),
828 ir->init_step*ir->delta_t);
832 fprintf(stderr, "%s steps, %s ps.\n",
833 gmx_step_str(ir->nsteps, sbuf), tbuf);
836 fprintf(fplog, "\n");
839 walltime_accounting_start(walltime_accounting);
840 wallcycle_start(wcycle, ewcRUN);
841 print_start(fplog, cr, walltime_accounting, "mdrun");
843 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
845 chkpt_ret = fcCheckPointParallel( cr->nodeid,
849 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
853 /***********************************************************
857 ************************************************************/
859 /* if rerunMD then read coordinates and velocities from input trajectory */
862 if (getenv("GMX_FORCE_UPDATE"))
870 bLastStep = !read_first_frame(oenv, &status,
871 opt2fn("-rerun", nfile, fnm),
872 &rerun_fr, TRX_NEED_X | TRX_READ_V);
873 if (rerun_fr.natoms != top_global->natoms)
876 "Number of atoms in trajectory (%d) does not match the "
877 "run input file (%d)\n",
878 rerun_fr.natoms, top_global->natoms);
880 if (ir->ePBC != epbcNONE)
884 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);
886 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
888 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
895 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
898 if (ir->ePBC != epbcNONE)
900 /* Set the shift vectors.
901 * Necessary here when have a static box different from the tpr box.
903 calc_shifts(rerun_fr.box, fr->shift_vec);
907 /* Loop over MD steps or if rerunMD to end of input trajectory,
908 * or, if max_hours>0, until max_hours is reached.
910 real max_hours = mdrunOptions.maximumHoursToRun;
912 /* Skip the first Nose-Hoover integration when we get the state from tpx */
913 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
914 bSumEkinhOld = FALSE;
916 bNeedRepartition = FALSE;
917 // TODO This implementation of ensemble orientation restraints is nasty because
918 // a user can't just do multi-sim with single-sim orientation restraints.
919 bUsingEnsembleRestraints = (fcd->disres.nsystems > 1) || (cr->ms && fcd->orires.nr);
922 // Replica exchange and ensemble restraints need all
923 // simulations to remain synchronized, so they need
924 // checkpoints and stop conditions to act on the same step, so
925 // the propagation of such signals must take place between
926 // simulations, not just within simulations.
927 bool checkpointIsLocal = !useReplicaExchange && !bUsingEnsembleRestraints;
928 bool stopConditionIsLocal = !useReplicaExchange && !bUsingEnsembleRestraints;
929 bool resetCountersIsLocal = true;
930 signals[eglsCHKPT] = SimulationSignal(checkpointIsLocal);
931 signals[eglsSTOPCOND] = SimulationSignal(stopConditionIsLocal);
932 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
935 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion = (DOMAINDECOMP(cr) ? DdOpenBalanceRegionBeforeForceComputation::yes : DdOpenBalanceRegionBeforeForceComputation::no);
936 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion = (DOMAINDECOMP(cr) ? DdCloseBalanceRegionAfterForceComputation::yes : DdCloseBalanceRegionAfterForceComputation::no);
938 step = ir->init_step;
941 // TODO extract this to new multi-simulation module
942 if (MASTER(cr) && MULTISIM(cr) && !useReplicaExchange)
944 if (!multisim_int_all_are_equal(cr->ms, ir->nsteps))
946 GMX_LOG(mdlog.warning).appendText(
947 "Note: The number of steps is not consistent across multi simulations,\n"
948 "but we are proceeding anyway!");
950 if (!multisim_int_all_are_equal(cr->ms, ir->init_step))
952 GMX_LOG(mdlog.warning).appendText(
953 "Note: The initial step is not consistent across multi simulations,\n"
954 "but we are proceeding anyway!");
958 /* and stop now if we should */
959 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
963 /* Determine if this is a neighbor search step */
964 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
966 if (bPMETune && bNStList)
968 /* PME grid + cut-off optimization with GPUs or PME nodes */
969 pme_loadbal_do(pme_loadbal, cr,
970 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
978 wallcycle_start(wcycle, ewcSTEP);
984 step = rerun_fr.step;
985 step_rel = step - ir->init_step;
998 bLastStep = (step_rel == ir->nsteps);
999 t = t0 + step*ir->delta_t;
1002 // TODO Refactor this, so that nstfep does not need a default value of zero
1003 if (ir->efep != efepNO || ir->bSimTemp)
1005 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
1006 requiring different logic. */
1011 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
1016 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
1018 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
1019 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
1020 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
1021 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
1024 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep &&
1025 do_per_step(step, replExParams.exchangeInterval));
1029 update_annealing_target_temp(ir, t, upd);
1032 if (bRerunMD && MASTER(cr))
1034 const bool constructVsites = (vsite && mdrunOptions.rerunConstructVsites);
1035 if (constructVsites && DOMAINDECOMP(cr))
1037 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
1039 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top->idef, ir->delta_t, *fr, graph, &bRerunWarnNoV);
1042 /* Stop Center of Mass motion */
1043 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
1047 /* for rerun MD always do Neighbour Searching */
1048 bNS = (bFirstStep || ir->nstlist != 0);
1052 /* Determine whether or not to do Neighbour Searching */
1053 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
1056 /* < 0 means stop at next step, > 0 means stop at next NS step */
1057 if ( (signals[eglsSTOPCOND].set < 0) ||
1058 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
1063 /* Determine whether or not to update the Born radii if doing GB */
1064 bBornRadii = bFirstStep;
1065 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
1070 /* do_log triggers energy and virial calculation. Because this leads
1071 * to different code paths, forces can be different. Thus for exact
1072 * continuation we should avoid extra log output.
1073 * Note that the || bLastStep can result in non-exact continuation
1074 * beyond the last step. But we don't consider that to be an issue.
1076 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
1077 do_verbose = mdrunOptions.verbose &&
1078 (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep || bRerunMD);
1080 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
1084 bMasterState = TRUE;
1088 bMasterState = FALSE;
1089 /* Correct the new box if it is too skewed */
1090 if (inputrecDynamicBox(ir))
1092 if (correct_box(fplog, step, state->box, graph))
1094 bMasterState = TRUE;
1097 if (DOMAINDECOMP(cr) && bMasterState)
1099 dd_collect_state(cr->dd, state, state_global);
1103 if (DOMAINDECOMP(cr))
1105 /* Repartition the domain decomposition */
1106 dd_partition_system(fplog, step, cr,
1107 bMasterState, nstglobalcomm,
1108 state_global, top_global, ir,
1109 state, &f, mdAtoms, top, fr,
1112 do_verbose && !bPMETunePrinting);
1113 shouldCheckNumberOfBondedInteractions = true;
1114 update_realloc(upd, state->natoms);
1118 if (MASTER(cr) && do_log)
1120 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
1123 if (ir->efep != efepNO)
1125 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1128 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1131 /* We need the kinetic energy at minus the half step for determining
1132 * the full step kinetic energy and possibly for T-coupling.*/
1133 /* This may not be quite working correctly yet . . . . */
1134 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1135 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1136 constr, &nullSignaller, state->box,
1137 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1138 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1139 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1140 top_global, top, state,
1141 &shouldCheckNumberOfBondedInteractions);
1143 clear_mat(force_vir);
1145 /* We write a checkpoint at this MD step when:
1146 * either at an NS step when we signalled through gs,
1147 * or at the last step (but not when we do not want confout),
1148 * but never at the first step or with rerun.
1150 bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1151 (bLastStep && mdrunOptions.writeConfout)) &&
1152 step > ir->init_step && !bRerunMD);
1155 signals[eglsCHKPT].set = 0;
1158 /* Determine the energy and pressure:
1159 * at nstcalcenergy steps and at energy output steps (set below).
1161 if (EI_VV(ir->eI) && (!bInitStep))
1163 /* for vv, the first half of the integration actually corresponds
1164 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1165 but the virial needs to be calculated on both the current step and the 'next' step. Future
1166 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1168 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1169 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1170 bCalcVir = bCalcEnerStep ||
1171 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1175 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1176 bCalcVir = bCalcEnerStep ||
1177 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1179 bCalcEner = bCalcEnerStep;
1181 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1183 if (do_ene || do_log || bDoReplEx)
1189 /* Do we need global communication ? */
1190 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1191 do_per_step(step, nstglobalcomm) ||
1192 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1194 force_flags = (GMX_FORCE_STATECHANGED |
1195 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1196 GMX_FORCE_ALLFORCES |
1197 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1198 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1199 (bDoFEP ? GMX_FORCE_DHDL : 0)
1204 /* Now is the time to relax the shells */
1205 relax_shell_flexcon(fplog, cr, mdrunOptions.verbose, step,
1206 ir, bNS, force_flags, top,
1208 state, &f, force_vir, mdatoms,
1209 nrnb, wcycle, graph, groups,
1210 shellfc, fr, bBornRadii, t, mu_tot,
1212 ddOpenBalanceRegion, ddCloseBalanceRegion);
1216 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
1217 (or the AWH update will be performed twice for one step when continuing). It would be best to
1218 call this update function from do_md_trajectory_writing but that would occur after do_force.
1219 One would have to divide the update_awh function into one function applying the AWH force
1220 and one doing the AWH bias update. The update AWH bias function could then be called after
1221 do_md_trajectory_writing (then containing update_awh_history).
1222 The checkpointing will in the future probably moved to the start of the md loop which will
1223 rid of this issue. */
1224 if (ir->bDoAwh && bCPT && MASTER(cr))
1226 ir->awh->updateHistory(state_global->awhHistory.get());
1229 /* The coordinates (x) are shifted (to get whole molecules)
1231 * This is parallellized as well, and does communication too.
1232 * Check comments in sim_util.c
1234 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1235 state->box, state->x, &state->hist,
1236 f, force_vir, mdatoms, enerd, fcd,
1237 state->lambda, graph,
1238 fr, vsite, mu_tot, t, ed, bBornRadii,
1239 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1240 ddOpenBalanceRegion, ddCloseBalanceRegion);
1243 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1244 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1246 rvec *vbuf = nullptr;
1248 wallcycle_start(wcycle, ewcUPDATE);
1249 if (ir->eI == eiVV && bInitStep)
1251 /* if using velocity verlet with full time step Ekin,
1252 * take the first half step only to compute the
1253 * virial for the first step. From there,
1254 * revert back to the initial coordinates
1255 * so that the input is actually the initial step.
1257 snew(vbuf, state->natoms);
1258 copy_rvecn(as_rvec_array(state->v.data()), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1262 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1263 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1266 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1267 ekind, M, upd, etrtVELOCITY1,
1270 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1272 wallcycle_stop(wcycle, ewcUPDATE);
1273 update_constraints(fplog, step, nullptr, ir, mdatoms,
1274 state, fr->bMolPBC, graph, f,
1275 &top->idef, shake_vir,
1276 cr, nrnb, wcycle, upd, constr,
1278 wallcycle_start(wcycle, ewcUPDATE);
1282 /* Need to unshift here if a do_force has been
1283 called in the previous step */
1284 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1286 /* if VV, compute the pressure and constraints */
1287 /* For VV2, we strictly only need this if using pressure
1288 * control, but we really would like to have accurate pressures
1290 * Think about ways around this in the future?
1291 * For now, keep this choice in comments.
1293 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1294 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1296 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1297 if (bCalcEner && ir->eI == eiVVAK)
1299 bSumEkinhOld = TRUE;
1301 /* for vv, the first half of the integration actually corresponds to the previous step.
1302 So we need information from the last step in the first half of the integration */
1303 if (bGStat || do_per_step(step-1, nstglobalcomm))
1305 wallcycle_stop(wcycle, ewcUPDATE);
1306 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1307 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1308 constr, &nullSignaller, state->box,
1309 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1310 (bGStat ? CGLO_GSTAT : 0)
1312 | (bTemp ? CGLO_TEMPERATURE : 0)
1313 | (bPres ? CGLO_PRESSURE : 0)
1314 | (bPres ? CGLO_CONSTRAINT : 0)
1315 | (bStopCM ? CGLO_STOPCM : 0)
1316 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1319 /* explanation of above:
1320 a) We compute Ekin at the full time step
1321 if 1) we are using the AveVel Ekin, and it's not the
1322 initial step, or 2) if we are using AveEkin, but need the full
1323 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1324 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1325 EkinAveVel because it's needed for the pressure */
1326 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1327 top_global, top, state,
1328 &shouldCheckNumberOfBondedInteractions);
1329 wallcycle_start(wcycle, ewcUPDATE);
1331 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1336 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1337 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1339 /* TODO This is only needed when we're about to write
1340 * a checkpoint, because we use it after the restart
1341 * (in a kludge?). But what should we be doing if
1342 * startingFromCheckpoint or bInitStep are true? */
1343 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1345 copy_mat(shake_vir, state->svir_prev);
1346 copy_mat(force_vir, state->fvir_prev);
1348 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1350 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1351 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1352 enerd->term[F_EKIN] = trace(ekind->ekin);
1355 else if (bExchanged)
1357 wallcycle_stop(wcycle, ewcUPDATE);
1358 /* We need the kinetic energy at minus the half step for determining
1359 * the full step kinetic energy and possibly for T-coupling.*/
1360 /* This may not be quite working correctly yet . . . . */
1361 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1362 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1363 constr, &nullSignaller, state->box,
1364 nullptr, &bSumEkinhOld,
1365 CGLO_GSTAT | CGLO_TEMPERATURE);
1366 wallcycle_start(wcycle, ewcUPDATE);
1369 /* if it's the initial step, we performed this first step just to get the constraint virial */
1370 if (ir->eI == eiVV && bInitStep)
1372 copy_rvecn(vbuf, as_rvec_array(state->v.data()), 0, state->natoms);
1375 wallcycle_stop(wcycle, ewcUPDATE);
1378 /* compute the conserved quantity */
1381 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1384 last_ekin = enerd->term[F_EKIN];
1386 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1388 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1390 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1391 if (ir->efep != efepNO && !bRerunMD)
1393 sum_dhdl(enerd, state->lambda, ir->fepvals);
1397 /* ######## END FIRST UPDATE STEP ############## */
1398 /* ######## If doing VV, we now have v(dt) ###### */
1401 /* perform extended ensemble sampling in lambda - we don't
1402 actually move to the new state before outputting
1403 statistics, but if performing simulated tempering, we
1404 do update the velocities and the tau_t. */
1406 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, as_rvec_array(state->v.data()), mdatoms);
1407 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1410 copy_df_history(state_global->dfhist, state->dfhist);
1414 /* Now we have the energies and forces corresponding to the
1415 * coordinates at time t. We must output all of this before
1418 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1419 ir, state, state_global, observablesHistory,
1421 outf, mdebin, ekind, f,
1423 bCPT, bRerunMD, bLastStep,
1424 mdrunOptions.writeConfout,
1426 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1427 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, as_rvec_array(state->x.data()), ir, t, wcycle);
1429 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1430 if (startingFromCheckpoint && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1432 copy_mat(state->svir_prev, shake_vir);
1433 copy_mat(state->fvir_prev, force_vir);
1436 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1438 /* Check whether everything is still allright */
1439 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1445 int nsteps_stop = -1;
1447 /* this just makes signals[].sig compatible with the hack
1448 of sending signals around by MPI_Reduce together with
1450 if ((gmx_get_stop_condition() == gmx_stop_cond_next_ns) ||
1451 (mdrunOptions.reproducible &&
1452 gmx_get_stop_condition() == gmx_stop_cond_next))
1454 /* We need at least two global communication steps to pass
1455 * around the signal. We stop at a pair-list creation step
1456 * to allow for exact continuation, when possible.
1458 signals[eglsSTOPCOND].sig = 1;
1459 nsteps_stop = std::max(ir->nstlist, 2*nstglobalcomm);
1461 else if (gmx_get_stop_condition() == gmx_stop_cond_next)
1463 /* Stop directly after the next global communication step.
1464 * This breaks exact continuation.
1466 signals[eglsSTOPCOND].sig = -1;
1467 nsteps_stop = nstglobalcomm + 1;
1472 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1473 gmx_get_signal_name(), nsteps_stop);
1477 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1478 gmx_get_signal_name(), nsteps_stop);
1480 handled_stop_condition = (int)gmx_get_stop_condition();
1482 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1483 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1484 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1486 /* Signal to terminate the run */
1487 signals[eglsSTOPCOND].sig = 1;
1490 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1492 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1495 if (bResetCountersHalfMaxH && MASTER(cr) &&
1496 elapsed_time > max_hours*60.0*60.0*0.495)
1498 /* Set flag that will communicate the signal to all ranks in the simulation */
1499 signals[eglsRESETCOUNTERS].sig = 1;
1502 /* In parallel we only have to check for checkpointing in steps
1503 * where we do global communication,
1504 * otherwise the other nodes don't know.
1506 const real cpt_period = mdrunOptions.checkpointOptions.period;
1507 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1510 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1511 signals[eglsCHKPT].set == 0)
1513 signals[eglsCHKPT].sig = 1;
1516 /* ######### START SECOND UPDATE STEP ################# */
1518 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1521 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1523 gmx_bool bIfRandomize;
1524 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1525 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1526 if (constr && bIfRandomize)
1528 update_constraints(fplog, step, nullptr, ir, mdatoms,
1529 state, fr->bMolPBC, graph, f,
1530 &top->idef, tmp_vir,
1531 cr, nrnb, wcycle, upd, constr,
1535 /* Box is changed in update() when we do pressure coupling,
1536 * but we should still use the old box for energy corrections and when
1537 * writing it to the energy file, so it matches the trajectory files for
1538 * the same timestep above. Make a copy in a separate array.
1540 copy_mat(state->box, lastbox);
1544 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1546 wallcycle_start(wcycle, ewcUPDATE);
1547 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1550 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1551 /* We can only do Berendsen coupling after we have summed
1552 * the kinetic energy or virial. Since the happens
1553 * in global_state after update, we should only do it at
1554 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1559 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1560 update_pcouple_before_coordinates(fplog, step, ir, state,
1561 parrinellorahmanMu, M,
1567 /* velocity half-step update */
1568 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1569 ekind, M, upd, etrtVELOCITY2,
1573 /* Above, initialize just copies ekinh into ekin,
1574 * it doesn't copy position (for VV),
1575 * and entire integrator for MD.
1578 if (ir->eI == eiVVAK)
1580 /* We probably only need md->homenr, not state->natoms */
1581 if (state->natoms > cbuf_nalloc)
1583 cbuf_nalloc = state->natoms;
1584 srenew(cbuf, cbuf_nalloc);
1586 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1589 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1590 ekind, M, upd, etrtPOSITION, cr, constr);
1591 wallcycle_stop(wcycle, ewcUPDATE);
1593 update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1594 fr->bMolPBC, graph, f,
1595 &top->idef, shake_vir,
1596 cr, nrnb, wcycle, upd, constr,
1599 if (ir->eI == eiVVAK)
1601 /* erase F_EKIN and F_TEMP here? */
1602 /* just compute the kinetic energy at the half step to perform a trotter step */
1603 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1604 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1605 constr, &nullSignaller, lastbox,
1606 nullptr, &bSumEkinhOld,
1607 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1609 wallcycle_start(wcycle, ewcUPDATE);
1610 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1611 /* now we know the scaling, we can compute the positions again again */
1612 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1614 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1615 ekind, M, upd, etrtPOSITION, cr, constr);
1616 wallcycle_stop(wcycle, ewcUPDATE);
1618 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1619 /* are the small terms in the shake_vir here due
1620 * to numerical errors, or are they important
1621 * physically? I'm thinking they are just errors, but not completely sure.
1622 * For now, will call without actually constraining, constr=NULL*/
1623 update_constraints(fplog, step, nullptr, ir, mdatoms,
1624 state, fr->bMolPBC, graph, f,
1625 &top->idef, tmp_vir,
1626 cr, nrnb, wcycle, upd, nullptr,
1631 /* this factor or 2 correction is necessary
1632 because half of the constraint force is removed
1633 in the vv step, so we have to double it. See
1634 the Redmine issue #1255. It is not yet clear
1635 if the factor of 2 is exact, or just a very
1636 good approximation, and this will be
1637 investigated. The next step is to see if this
1638 can be done adding a dhdl contribution from the
1639 rattle step, but this is somewhat more
1640 complicated with the current code. Will be
1641 investigated, hopefully for 4.6.3. However,
1642 this current solution is much better than
1643 having it completely wrong.
1645 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1649 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1654 /* Need to unshift here */
1655 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1658 if (vsite != nullptr)
1660 wallcycle_start(wcycle, ewcVSITECONSTR);
1661 if (graph != nullptr)
1663 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1665 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
1666 top->idef.iparams, top->idef.il,
1667 fr->ePBC, fr->bMolPBC, cr, state->box);
1669 if (graph != nullptr)
1671 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1673 wallcycle_stop(wcycle, ewcVSITECONSTR);
1676 /* ############## IF NOT VV, Calculate globals HERE ############ */
1677 /* With Leap-Frog we can skip compute_globals at
1678 * non-communication steps, but we need to calculate
1679 * the kinetic energy one step before communication.
1682 // Organize to do inter-simulation signalling on steps if
1683 // and when algorithms require it.
1684 bool doInterSimSignal = (!bFirstStep && bDoReplEx) || bUsingEnsembleRestraints;
1686 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1688 // Since we're already communicating at this step, we
1689 // can propagate intra-simulation signals. Note that
1690 // check_nstglobalcomm has the responsibility for
1691 // choosing the value of nstglobalcomm that is one way
1692 // bGStat becomes true, so we can't get into a
1693 // situation where e.g. checkpointing can't be
1695 bool doIntraSimSignal = true;
1696 SimulationSignaller signaller(&signals, cr, doInterSimSignal, doIntraSimSignal);
1698 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1699 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1702 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1703 (bGStat ? CGLO_GSTAT : 0)
1704 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1705 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1706 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1707 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1709 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1711 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1712 top_global, top, state,
1713 &shouldCheckNumberOfBondedInteractions);
1717 /* ############# END CALC EKIN AND PRESSURE ################# */
1719 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1720 the virial that should probably be addressed eventually. state->veta has better properies,
1721 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1722 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1724 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1726 /* Sum up the foreign energy and dhdl terms for md and sd.
1727 Currently done every step so that dhdl is correct in the .edr */
1728 sum_dhdl(enerd, state->lambda, ir->fepvals);
1731 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1732 pres, force_vir, shake_vir,
1736 /* ################# END UPDATE STEP 2 ################# */
1737 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1739 /* The coordinates (x) were unshifted in update */
1742 /* We will not sum ekinh_old,
1743 * so signal that we still have to do it.
1745 bSumEkinhOld = TRUE;
1750 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1752 /* use the directly determined last velocity, not actually the averaged half steps */
1753 if (bTrotter && ir->eI == eiVV)
1755 enerd->term[F_EKIN] = last_ekin;
1757 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1759 if (integratorHasConservedEnergyQuantity(ir))
1763 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1767 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1770 /* ######### END PREPARING EDR OUTPUT ########### */
1776 if (fplog && do_log && bDoExpanded)
1778 /* only needed if doing expanded ensemble */
1779 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1780 state_global->dfhist, state->fep_state, ir->nstlog, step);
1784 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1785 t, mdatoms->tmass, enerd, state,
1786 ir->fepvals, ir->expandedvals, lastbox,
1787 shake_vir, force_vir, total_vir, pres,
1788 ekind, mu_tot, constr);
1792 upd_mdebin_step(mdebin);
1795 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1796 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1798 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1800 eprNORMAL, mdebin, fcd, groups, &(ir->opts), ir->awh);
1804 pull_print_output(ir->pull_work, step, t);
1807 if (do_per_step(step, ir->nstlog))
1809 if (fflush(fplog) != 0)
1811 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1817 /* Have to do this part _after_ outputting the logfile and the edr file */
1818 /* Gets written into the state at the beginning of next loop*/
1819 state->fep_state = lamnew;
1821 /* Print the remaining wall clock time for the run */
1822 if (MULTIMASTER(cr) &&
1823 (do_verbose || gmx_got_usr_signal()) &&
1828 fprintf(stderr, "\n");
1830 print_time(stderr, walltime_accounting, step, ir, cr);
1833 /* Ion/water position swapping.
1834 * Not done in last step since trajectory writing happens before this call
1835 * in the MD loop and exchanges would be lost anyway. */
1836 bNeedRepartition = FALSE;
1837 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1838 do_per_step(step, ir->swap->nstswap))
1840 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1841 bRerunMD ? rerun_fr.x : as_rvec_array(state->x.data()),
1842 bRerunMD ? rerun_fr.box : state->box,
1843 MASTER(cr) && mdrunOptions.verbose,
1846 if (bNeedRepartition && DOMAINDECOMP(cr))
1848 dd_collect_state(cr->dd, state, state_global);
1852 /* Replica exchange */
1856 bExchanged = replica_exchange(fplog, cr, repl_ex,
1857 state_global, enerd,
1861 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1863 dd_partition_system(fplog, step, cr, TRUE, 1,
1864 state_global, top_global, ir,
1865 state, &f, mdAtoms, top, fr,
1867 nrnb, wcycle, FALSE);
1868 shouldCheckNumberOfBondedInteractions = true;
1869 update_realloc(upd, state->natoms);
1874 startingFromCheckpoint = false;
1876 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1877 /* With all integrators, except VV, we need to retain the pressure
1878 * at the current step for coupling at the next step.
1880 if ((state->flags & (1<<estPRES_PREV)) &&
1882 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1884 /* Store the pressure in t_state for pressure coupling
1885 * at the next MD step.
1887 copy_mat(pres, state->pres_prev);
1890 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1892 if ( (membed != nullptr) && (!bLastStep) )
1894 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1901 /* read next frame from input trajectory */
1902 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1907 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1911 cycles = wallcycle_stop(wcycle, ewcSTEP);
1912 if (DOMAINDECOMP(cr) && wcycle)
1914 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1917 if (!bRerunMD || !rerun_fr.bStep)
1919 /* increase the MD step number */
1924 /* TODO make a counter-reset module */
1925 /* If it is time to reset counters, set a flag that remains
1926 true until counters actually get reset */
1927 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1928 signals[eglsRESETCOUNTERS].set != 0)
1930 if (pme_loadbal_is_active(pme_loadbal))
1932 /* Do not permit counter reset while PME load
1933 * balancing is active. The only purpose for resetting
1934 * counters is to measure reliable performance data,
1935 * and that can't be done before balancing
1938 * TODO consider fixing this by delaying the reset
1939 * until after load balancing completes,
1940 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1941 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1942 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1943 "resetting counters later in the run, e.g. with gmx "
1944 "mdrun -resetstep.", step);
1946 reset_all_counters(fplog, mdlog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1947 use_GPU(fr->nbv) ? fr->nbv : nullptr, fr->pmedata);
1948 wcycle_set_reset_counters(wcycle, -1);
1949 if (!thisRankHasDuty(cr, DUTY_PME))
1951 /* Tell our PME node to reset its counters */
1952 gmx_pme_send_resetcounters(cr, step);
1954 /* Correct max_hours for the elapsed time */
1955 max_hours -= elapsed_time/(60.0*60.0);
1956 /* If mdrun -maxh -resethway was active, it can only trigger once */
1957 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1958 /* Reset can only happen once, so clear the triggering flag. */
1959 signals[eglsRESETCOUNTERS].set = 0;
1962 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1963 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1966 /* End of main MD loop */
1968 /* Closing TNG files can include compressing data. Therefore it is good to do that
1969 * before stopping the time measurements. */
1970 mdoutf_tng_close(outf);
1972 /* Stop measuring walltime */
1973 walltime_accounting_end(walltime_accounting);
1975 if (bRerunMD && MASTER(cr))
1980 if (!thisRankHasDuty(cr, DUTY_PME))
1982 /* Tell the PME only node to finish */
1983 gmx_pme_send_finish(cr);
1988 if (ir->nstcalcenergy > 0 && !bRerunMD)
1990 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1991 eprAVER, mdebin, fcd, groups, &(ir->opts), ir->awh);
1999 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
2002 done_shellfc(fplog, shellfc, step_rel);
2004 if (useReplicaExchange && MASTER(cr))
2006 print_replica_exchange_statistics(fplog, repl_ex);
2014 // Clean up swapcoords
2015 if (ir->eSwapCoords != eswapNO)
2017 finish_swapcoords(ir->swap);
2020 /* Do essential dynamics cleanup if needed. Close .edo file */
2023 /* IMD cleanup, if bIMD is TRUE. */
2024 IMD_finalize(ir->bIMD, ir->imd);
2026 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
2027 if (step_rel >= wcycle_get_reset_counters(wcycle) &&
2028 signals[eglsRESETCOUNTERS].set == 0 &&
2029 !bResetCountersHalfMaxH)
2031 walltime_accounting_set_valid_finish(walltime_accounting);