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.
39 * \brief Implements the integrator for normal molecular dynamics simulations
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
56 #include "gromacs/awh/awh.h"
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/compat/make_unique.h"
59 #include "gromacs/domdec/collect.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_network.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/essentialdynamics/edsam.h"
64 #include "gromacs/ewald/pme.h"
65 #include "gromacs/ewald/pme-load-balancing.h"
66 #include "gromacs/fileio/trxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/gpu_utils/gpu_utils.h"
70 #include "gromacs/imd/imd.h"
71 #include "gromacs/listed-forces/manage-threading.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/utilities.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/math/vectypes.h"
76 #include "gromacs/mdlib/compute_io.h"
77 #include "gromacs/mdlib/constr.h"
78 #include "gromacs/mdlib/ebin.h"
79 #include "gromacs/mdlib/expanded.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/md_support.h"
84 #include "gromacs/mdlib/mdatoms.h"
85 #include "gromacs/mdlib/mdebin.h"
86 #include "gromacs/mdlib/mdoutf.h"
87 #include "gromacs/mdlib/mdrun.h"
88 #include "gromacs/mdlib/mdsetup.h"
89 #include "gromacs/mdlib/membed.h"
90 #include "gromacs/mdlib/nb_verlet.h"
91 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
92 #include "gromacs/mdlib/ns.h"
93 #include "gromacs/mdlib/repl_ex.h"
94 #include "gromacs/mdlib/shellfc.h"
95 #include "gromacs/mdlib/sighandler.h"
96 #include "gromacs/mdlib/sim_util.h"
97 #include "gromacs/mdlib/simulationsignal.h"
98 #include "gromacs/mdlib/tgroup.h"
99 #include "gromacs/mdlib/trajectory_writing.h"
100 #include "gromacs/mdlib/update.h"
101 #include "gromacs/mdlib/vcm.h"
102 #include "gromacs/mdlib/vsite.h"
103 #include "gromacs/mdtypes/awh-history.h"
104 #include "gromacs/mdtypes/awh-params.h"
105 #include "gromacs/mdtypes/commrec.h"
106 #include "gromacs/mdtypes/df_history.h"
107 #include "gromacs/mdtypes/energyhistory.h"
108 #include "gromacs/mdtypes/fcdata.h"
109 #include "gromacs/mdtypes/forcerec.h"
110 #include "gromacs/mdtypes/group.h"
111 #include "gromacs/mdtypes/inputrec.h"
112 #include "gromacs/mdtypes/interaction_const.h"
113 #include "gromacs/mdtypes/md_enums.h"
114 #include "gromacs/mdtypes/mdatom.h"
115 #include "gromacs/mdtypes/observableshistory.h"
116 #include "gromacs/mdtypes/state.h"
117 #include "gromacs/pbcutil/mshift.h"
118 #include "gromacs/pbcutil/pbc.h"
119 #include "gromacs/pulling/pull.h"
120 #include "gromacs/swap/swapcoords.h"
121 #include "gromacs/timing/wallcycle.h"
122 #include "gromacs/timing/walltime_accounting.h"
123 #include "gromacs/topology/atoms.h"
124 #include "gromacs/topology/idef.h"
125 #include "gromacs/topology/mtop_util.h"
126 #include "gromacs/topology/topology.h"
127 #include "gromacs/trajectory/trajectoryframe.h"
128 #include "gromacs/utility/basedefinitions.h"
129 #include "gromacs/utility/cstringutil.h"
130 #include "gromacs/utility/fatalerror.h"
131 #include "gromacs/utility/logger.h"
132 #include "gromacs/utility/real.h"
133 #include "gromacs/utility/smalloc.h"
135 #include "integrator.h"
138 #include "corewrap.h"
141 using gmx::SimulationSignaller;
143 //! Resets all the counters.
144 static void reset_all_counters(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
146 gmx_int64_t *step_rel, t_inputrec *ir,
147 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
148 gmx_walltime_accounting_t walltime_accounting,
149 struct nonbonded_verlet_t *nbv,
150 struct gmx_pme_t *pme)
152 char sbuf[STEPSTRSIZE];
154 /* Reset all the counters related to performance over the run */
155 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
156 "step %s: resetting all time and cycle counters",
157 gmx_step_str(step, sbuf));
161 nbnxn_gpu_reset_timings(nbv);
164 if (pme_gpu_task_enabled(pme))
166 pme_gpu_reset_timings(pme);
169 if (use_GPU(nbv) || pme_gpu_task_enabled(pme))
174 wallcycle_stop(wcycle, ewcRUN);
175 wallcycle_reset_all(wcycle);
176 if (DOMAINDECOMP(cr))
178 reset_dd_statistics_counters(cr->dd);
181 ir->init_step += *step_rel;
182 ir->nsteps -= *step_rel;
184 wallcycle_start(wcycle, ewcRUN);
185 walltime_accounting_start(walltime_accounting);
186 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
189 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
191 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
192 * \param[in,out] globalState The global state container
193 * \param[in] constructVsites When true, vsite coordinates are constructed
194 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
195 * \param[in] idef Topology parameters, used for constructing vsites
196 * \param[in] timeStep Time step, used for constructing vsites
197 * \param[in] forceRec Force record, used for constructing vsites
198 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
199 * \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
201 static void prepareRerunState(const t_trxframe &rerunFrame,
202 t_state *globalState,
203 bool constructVsites,
204 const gmx_vsite_t *vsite,
207 const t_forcerec &forceRec,
209 gmx_bool *warnWhenNoV)
211 for (int i = 0; i < globalState->natoms; i++)
213 copy_rvec(rerunFrame.x[i], globalState->x[i]);
217 for (int i = 0; i < globalState->natoms; i++)
219 copy_rvec(rerunFrame.v[i], globalState->v[i]);
224 for (int i = 0; i < globalState->natoms; i++)
226 clear_rvec(globalState->v[i]);
230 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
231 " Ekin, temperature and pressure are incorrect,\n"
232 " the virial will be incorrect when constraints are present.\n"
234 *warnWhenNoV = FALSE;
237 copy_mat(rerunFrame.box, globalState->box);
241 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
245 /* Following is necessary because the graph may get out of sync
246 * with the coordinates if we only have every N'th coordinate set
248 mk_mshift(nullptr, graph, forceRec.ePBC, globalState->box, as_rvec_array(globalState->x.data()));
249 shift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
251 construct_vsites(vsite, as_rvec_array(globalState->x.data()), timeStep, as_rvec_array(globalState->v.data()),
252 idef.iparams, idef.il,
253 forceRec.ePBC, forceRec.bMolPBC, nullptr, globalState->box);
256 unshift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
261 void gmx::Integrator::do_md()
263 // TODO Historically, the EM and MD "integrators" used different
264 // names for the t_inputrec *parameter, but these must have the
265 // same name, now that it's a member of a struct. We use this ir
266 // alias to avoid a large ripple of nearly useless changes.
267 // t_inputrec is being replaced by IMdpOptionsProvider, so this
268 // will go away eventually.
269 t_inputrec *ir = inputrec;
270 gmx_mdoutf *outf = nullptr;
271 gmx_int64_t step, step_rel;
273 double t, t0, lam0[efptNR];
274 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
275 gmx_bool bNS, bNStList, bSimAnn, bStopCM,
276 bFirstStep, bInitStep, bLastStep = FALSE;
277 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
278 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
279 bForceUpdate = FALSE, bCPT;
280 gmx_bool bMasterState;
281 int force_flags, cglo_flags;
282 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
287 matrix parrinellorahmanMu, M;
289 gmx_repl_ex_t repl_ex = nullptr;
292 t_mdebin *mdebin = nullptr;
293 gmx_enerdata_t *enerd;
294 PaddedRVecVector f {};
295 gmx_global_stat_t gstat;
296 gmx_update_t *upd = nullptr;
297 t_graph *graph = nullptr;
298 gmx_groups_t *groups;
299 gmx_ekindata_t *ekind;
300 gmx_shellfc_t *shellfc;
301 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
302 gmx_bool bResetCountersHalfMaxH = FALSE;
303 gmx_bool bTemp, bPres, bTrotter;
305 rvec *cbuf = nullptr;
312 real saved_conserved_quantity = 0;
316 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
317 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
320 /* PME load balancing data for GPU kernels */
321 pme_load_balancing_t *pme_loadbal = nullptr;
322 gmx_bool bPMETune = FALSE;
323 gmx_bool bPMETunePrinting = FALSE;
326 gmx_bool bIMDstep = FALSE;
329 /* Temporary addition for FAHCORE checkpointing */
332 /* Domain decomposition could incorrectly miss a bonded
333 interaction, but checking for that requires a global
334 communication stage, which does not otherwise happen in DD
335 code. So we do that alongside the first global energy reduction
336 after a new DD is made. These variables handle whether the
337 check happens, and the result it returns. */
338 bool shouldCheckNumberOfBondedInteractions = false;
339 int totalNumberOfBondedInteractions = -1;
341 SimulationSignals signals;
342 // Most global communnication stages don't propagate mdrun
343 // signals, and will use this object to achieve that.
344 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
346 if (!mdrunOptions.writeConfout)
348 // This is on by default, and the main known use case for
349 // turning it off is for convenience in benchmarking, which is
350 // something that should not show up in the general user
352 GMX_LOG(mdlog.info).asParagraph().
353 appendText("The -noconfout functionality is deprecated, and may be removed in a future version.");
355 if (mdrunOptions.timingOptions.resetHalfway)
357 GMX_LOG(mdlog.info).asParagraph().
358 appendText("The -resethway functionality is deprecated, and may be removed in a future version.");
361 /* Signal to reset the counters half the simulation steps. */
362 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
364 /* Signal to reset the counters halfway the simulation time. */
365 bResetCountersHalfMaxH = (mdrunOptions.maximumHoursToRun > 0);
368 /* md-vv uses averaged full step velocities for T-control
369 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
370 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
371 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
373 const gmx_bool bRerunMD = mdrunOptions.rerun;
374 int nstglobalcomm = mdrunOptions.globalCommunicationInterval;
378 /* Since we don't know if the frames read are related in any way,
379 * rebuild the neighborlist at every step.
382 ir->nstcalcenergy = 1;
386 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
387 bGStatEveryStep = (nstglobalcomm == 1);
391 ir->nstxout_compressed = 0;
393 groups = &top_global->groups;
395 gmx_edsam *ed = nullptr;
396 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
398 /* Initialize essential dynamics sampling */
399 ed = init_edsam(opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm),
402 state_global, observablesHistory,
403 oenv, mdrunOptions.continuationOptions.appendFiles);
407 init_md(fplog, cr, outputProvider, ir, oenv, mdrunOptions,
408 &t, &t0, state_global, lam0,
409 nrnb, top_global, &upd, deform,
410 nfile, fnm, &outf, &mdebin,
411 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, wcycle);
413 clear_mat(total_vir);
415 /* Energy terms and groups */
417 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
420 /* Kinetic energy data */
422 init_ekindata(fplog, top_global, &(ir->opts), ekind);
423 /* Copy the cos acceleration to the groups struct */
424 ekind->cosacc.cos_accel = ir->cos_accel;
426 gstat = global_stat_init(ir);
428 /* Check for polarizable models and flexible constraints */
429 shellfc = init_shell_flexcon(fplog,
430 top_global, constr ? constr->numFlexibleConstraints() : 0,
431 ir->nstcalcenergy, DOMAINDECOMP(cr));
434 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
435 if ((io > 2000) && MASTER(cr))
438 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
443 /* Set up interactive MD (IMD) */
444 init_IMD(ir, cr, ms, top_global, fplog, ir->nstcalcenergy,
445 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
446 nfile, fnm, oenv, mdrunOptions);
448 // Local state only becomes valid now.
449 std::unique_ptr<t_state> stateInstance;
452 if (DOMAINDECOMP(cr))
454 top = dd_init_local_top(top_global);
456 stateInstance = compat::make_unique<t_state>();
457 state = stateInstance.get();
458 dd_init_local_state(cr->dd, state_global, state);
460 /* Distribute the charge groups over the nodes from the master node */
461 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
462 state_global, top_global, ir,
464 state, &f, mdAtoms, top, fr,
466 nrnb, nullptr, FALSE);
467 shouldCheckNumberOfBondedInteractions = true;
468 update_realloc(upd, state->natoms);
472 state_change_natoms(state_global, state_global->natoms);
473 /* We need to allocate one element extra, since we might use
474 * (unaligned) 4-wide SIMD loads to access rvec entries.
476 f.resize(gmx::paddedRVecVectorSize(state_global->natoms));
477 /* Copy the pointer to the global state */
478 state = state_global;
481 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
482 &graph, mdAtoms, constr, vsite, shellfc);
484 update_realloc(upd, state->natoms);
487 auto mdatoms = mdAtoms->mdatoms();
489 // NOTE: The global state is no longer used at this point.
490 // But state_global is still used as temporary storage space for writing
491 // the global state to file and potentially for replica exchange.
492 // (Global topology should persist.)
494 update_mdatoms(mdatoms, state->lambda[efptMASS]);
496 const ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
497 bool startingFromCheckpoint = continuationOptions.startedFromCheckpoint;
501 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
506 if (startingFromCheckpoint)
508 /* Update mdebin with energy history if appending to output files */
509 if (continuationOptions.appendFiles)
511 restore_energyhistory_from_state(mdebin, observablesHistory->energyHistory.get());
513 else if (observablesHistory->energyHistory != nullptr)
515 /* We might have read an energy history from checkpoint.
516 * As we are not appending, we want to restart the statistics.
517 * Free the allocated memory and reset the counts.
519 observablesHistory->energyHistory = {};
522 if (observablesHistory->energyHistory == nullptr)
524 observablesHistory->energyHistory = compat::make_unique<energyhistory_t>();
526 /* Set the initial energy history in state by updating once */
527 update_energyhistory(observablesHistory->energyHistory.get(), mdebin);
530 // TODO: Remove this by converting AWH into a ForceProvider
531 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms, startingFromCheckpoint,
533 opt2fn("-awh", nfile, fnm), ir->pull_work);
535 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
536 if (useReplicaExchange && MASTER(cr))
538 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir,
541 /* PME tuning is only supported in the Verlet scheme, with PME for
542 * Coulomb. It is not supported with only LJ PME, or for
544 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !bRerunMD &&
545 !mdrunOptions.reproducible && ir->cutoff_scheme != ecutsGROUP);
548 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
549 fr->ic, fr->nbv->listParams.get(), fr->pmedata, use_GPU(fr->nbv),
553 if (!ir->bContinuation && !bRerunMD)
555 if (state->flags & (1 << estV))
557 /* Set the velocities of vsites, shells and frozen atoms to zero */
558 for (i = 0; i < mdatoms->homenr; i++)
560 if (mdatoms->ptype[i] == eptVSite ||
561 mdatoms->ptype[i] == eptShell)
563 clear_rvec(state->v[i]);
565 else if (mdatoms->cFREEZE)
567 for (m = 0; m < DIM; m++)
569 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
580 /* Constrain the initial coordinates and velocities */
581 do_constrain_first(fplog, constr, ir, mdatoms, state);
585 /* Construct the virtual sites for the initial configuration */
586 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, nullptr,
587 top->idef.iparams, top->idef.il,
588 fr->ePBC, fr->bMolPBC, cr, state->box);
592 if (ir->efep != efepNO)
594 /* Set free energy calculation frequency as the greatest common
595 * denominator of nstdhdl and repl_ex_nst. */
596 nstfep = ir->fepvals->nstdhdl;
599 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
601 if (useReplicaExchange)
603 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
607 /* Be REALLY careful about what flags you set here. You CANNOT assume
608 * this is the first step, since we might be restarting from a checkpoint,
609 * and in that case we should not do any modifications to the state.
611 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
613 if (continuationOptions.haveReadEkin)
615 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
618 cglo_flags = (CGLO_INITIALIZATION | CGLO_TEMPERATURE | CGLO_GSTAT
619 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
620 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
621 | (continuationOptions.haveReadEkin ? CGLO_READEKIN : 0));
623 bSumEkinhOld = FALSE;
624 /* To minimize communication, compute_globals computes the COM velocity
625 * and the kinetic energy for the velocities without COM motion removed.
626 * Thus to get the kinetic energy without the COM contribution, we need
627 * to call compute_globals twice.
629 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
631 int cglo_flags_iteration = cglo_flags;
632 if (bStopCM && cgloIteration == 0)
634 cglo_flags_iteration |= CGLO_STOPCM;
635 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
637 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
638 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
639 constr, &nullSignaller, state->box,
640 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags_iteration
641 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
643 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
644 top_global, top, state,
645 &shouldCheckNumberOfBondedInteractions);
646 if (ir->eI == eiVVAK)
648 /* a second call to get the half step temperature initialized as well */
649 /* we do the same call as above, but turn the pressure off -- internally to
650 compute_globals, this is recognized as a velocity verlet half-step
651 kinetic energy calculation. This minimized excess variables, but
652 perhaps loses some logic?*/
654 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
655 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
656 constr, &nullSignaller, state->box,
657 nullptr, &bSumEkinhOld,
658 cglo_flags & ~CGLO_PRESSURE);
661 /* Calculate the initial half step temperature, and save the ekinh_old */
662 if (!continuationOptions.startedFromCheckpoint)
664 for (i = 0; (i < ir->opts.ngtc); i++)
666 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
670 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
671 temperature control */
672 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
676 if (!ir->bContinuation)
678 if (constr && ir->eConstrAlg == econtLINCS)
681 "RMS relative constraint deviation after constraining: %.2e\n",
684 if (EI_STATE_VELOCITY(ir->eI))
686 real temp = enerd->term[F_TEMP];
689 /* Result of Ekin averaged over velocities of -half
690 * and +half step, while we only have -half step here.
694 fprintf(fplog, "Initial temperature: %g K\n", temp);
700 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
701 " input trajectory '%s'\n\n",
702 *(top_global->name), opt2fn("-rerun", nfile, fnm));
703 if (mdrunOptions.verbose)
705 fprintf(stderr, "Calculated time to finish depends on nsteps from "
706 "run input file,\nwhich may not correspond to the time "
707 "needed to process input trajectory.\n\n");
713 fprintf(stderr, "starting mdrun '%s'\n",
714 *(top_global->name));
717 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
721 sprintf(tbuf, "%s", "infinite");
723 if (ir->init_step > 0)
725 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
726 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
727 gmx_step_str(ir->init_step, sbuf2),
728 ir->init_step*ir->delta_t);
732 fprintf(stderr, "%s steps, %s ps.\n",
733 gmx_step_str(ir->nsteps, sbuf), tbuf);
736 fprintf(fplog, "\n");
739 walltime_accounting_start(walltime_accounting);
740 wallcycle_start(wcycle, ewcRUN);
741 print_start(fplog, cr, walltime_accounting, "mdrun");
743 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
745 chkpt_ret = fcCheckPointParallel( cr->nodeid,
749 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
753 /***********************************************************
757 ************************************************************/
759 /* if rerunMD then read coordinates and velocities from input trajectory */
762 if (getenv("GMX_FORCE_UPDATE"))
770 bLastStep = !read_first_frame(oenv, &status,
771 opt2fn("-rerun", nfile, fnm),
772 &rerun_fr, TRX_NEED_X | TRX_READ_V);
773 if (rerun_fr.natoms != top_global->natoms)
776 "Number of atoms in trajectory (%d) does not match the "
777 "run input file (%d)\n",
778 rerun_fr.natoms, top_global->natoms);
780 if (ir->ePBC != epbcNONE)
784 gmx_fatal(FARGS, "Rerun trajectory frame step %ld time %f does not contain a box, while pbc is used", rerun_fr.step, rerun_fr.time);
786 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
788 gmx_fatal(FARGS, "Rerun trajectory frame step %ld time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
795 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
798 if (ir->ePBC != epbcNONE)
800 /* Set the shift vectors.
801 * Necessary here when have a static box different from the tpr box.
803 calc_shifts(rerun_fr.box, fr->shift_vec);
807 /* Loop over MD steps or if rerunMD to end of input trajectory,
808 * or, if max_hours>0, until max_hours is reached.
810 real max_hours = mdrunOptions.maximumHoursToRun;
812 /* Skip the first Nose-Hoover integration when we get the state from tpx */
813 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
814 bSumEkinhOld = FALSE;
816 bNeedRepartition = FALSE;
818 bool simulationsShareState = false;
819 int nstSignalComm = nstglobalcomm;
821 // TODO This implementation of ensemble orientation restraints is nasty because
822 // a user can't just do multi-sim with single-sim orientation restraints.
823 bool usingEnsembleRestraints = (fcd->disres.nsystems > 1) || (ms && fcd->orires.nr);
824 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && ms);
826 // Replica exchange, ensemble restraints and AWH need all
827 // simulations to remain synchronized, so they need
828 // checkpoints and stop conditions to act on the same step, so
829 // the propagation of such signals must take place between
830 // simulations, not just within simulations.
831 // TODO: Make algorithm initializers set these flags.
832 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
833 bool resetCountersIsLocal = true;
834 signals[eglsCHKPT] = SimulationSignal(!simulationsShareState);
835 signals[eglsSTOPCOND] = SimulationSignal(!simulationsShareState);
836 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
838 if (simulationsShareState)
840 // Inter-simulation signal communication does not need to happen
841 // often, so we use a minimum of 200 steps to reduce overhead.
842 const int c_minimumInterSimulationSignallingInterval = 200;
843 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1)/nstglobalcomm)*nstglobalcomm;
847 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion = (DOMAINDECOMP(cr) ? DdOpenBalanceRegionBeforeForceComputation::yes : DdOpenBalanceRegionBeforeForceComputation::no);
848 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion = (DOMAINDECOMP(cr) ? DdCloseBalanceRegionAfterForceComputation::yes : DdCloseBalanceRegionAfterForceComputation::no);
850 step = ir->init_step;
853 // TODO extract this to new multi-simulation module
854 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
856 if (!multisim_int_all_are_equal(ms, ir->nsteps))
858 GMX_LOG(mdlog.warning).appendText(
859 "Note: The number of steps is not consistent across multi simulations,\n"
860 "but we are proceeding anyway!");
862 if (!multisim_int_all_are_equal(ms, ir->init_step))
864 GMX_LOG(mdlog.warning).appendText(
865 "Note: The initial step is not consistent across multi simulations,\n"
866 "but we are proceeding anyway!");
870 /* and stop now if we should */
871 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
875 /* Determine if this is a neighbor search step */
876 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
878 if (bPMETune && bNStList)
880 /* PME grid + cut-off optimization with GPUs or PME nodes */
881 pme_loadbal_do(pme_loadbal, cr,
882 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
890 wallcycle_start(wcycle, ewcSTEP);
896 step = rerun_fr.step;
897 step_rel = step - ir->init_step;
910 bLastStep = (step_rel == ir->nsteps);
911 t = t0 + step*ir->delta_t;
914 // TODO Refactor this, so that nstfep does not need a default value of zero
915 if (ir->efep != efepNO || ir->bSimTemp)
917 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
918 requiring different logic. */
923 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
928 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
930 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
931 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
932 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
933 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
936 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep &&
937 do_per_step(step, replExParams.exchangeInterval));
941 update_annealing_target_temp(ir, t, upd);
944 if (bRerunMD && MASTER(cr))
946 const bool constructVsites = (vsite && mdrunOptions.rerunConstructVsites);
947 if (constructVsites && DOMAINDECOMP(cr))
949 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
951 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top->idef, ir->delta_t, *fr, graph, &bRerunWarnNoV);
954 /* Stop Center of Mass motion */
955 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
959 /* for rerun MD always do Neighbour Searching */
960 bNS = (bFirstStep || ir->nstlist != 0);
964 /* Determine whether or not to do Neighbour Searching */
965 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
968 /* < 0 means stop at next step, > 0 means stop at next NS step */
969 if ( (signals[eglsSTOPCOND].set < 0) ||
970 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
975 /* do_log triggers energy and virial calculation. Because this leads
976 * to different code paths, forces can be different. Thus for exact
977 * continuation we should avoid extra log output.
978 * Note that the || bLastStep can result in non-exact continuation
979 * beyond the last step. But we don't consider that to be an issue.
981 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
982 do_verbose = mdrunOptions.verbose &&
983 (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep || bRerunMD);
985 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
993 bMasterState = FALSE;
994 /* Correct the new box if it is too skewed */
995 if (inputrecDynamicBox(ir))
997 if (correct_box(fplog, step, state->box, graph))
1002 if (DOMAINDECOMP(cr) && bMasterState)
1004 dd_collect_state(cr->dd, state, state_global);
1008 if (DOMAINDECOMP(cr))
1010 /* Repartition the domain decomposition */
1011 dd_partition_system(fplog, step, cr,
1012 bMasterState, nstglobalcomm,
1013 state_global, top_global, ir,
1015 state, &f, mdAtoms, top, fr,
1018 do_verbose && !bPMETunePrinting);
1019 shouldCheckNumberOfBondedInteractions = true;
1020 update_realloc(upd, state->natoms);
1024 if (MASTER(cr) && do_log)
1026 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
1029 if (ir->efep != efepNO)
1031 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1034 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1037 /* We need the kinetic energy at minus the half step for determining
1038 * the full step kinetic energy and possibly for T-coupling.*/
1039 /* This may not be quite working correctly yet . . . . */
1040 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1041 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1042 constr, &nullSignaller, state->box,
1043 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1044 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1045 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1046 top_global, top, state,
1047 &shouldCheckNumberOfBondedInteractions);
1049 clear_mat(force_vir);
1051 /* We write a checkpoint at this MD step when:
1052 * either at an NS step when we signalled through gs,
1053 * or at the last step (but not when we do not want confout),
1054 * but never at the first step or with rerun.
1056 bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1057 (bLastStep && mdrunOptions.writeConfout)) &&
1058 step > ir->init_step && !bRerunMD);
1061 signals[eglsCHKPT].set = 0;
1064 /* Determine the energy and pressure:
1065 * at nstcalcenergy steps and at energy output steps (set below).
1067 if (EI_VV(ir->eI) && (!bInitStep))
1069 /* for vv, the first half of the integration actually corresponds
1070 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1071 but the virial needs to be calculated on both the current step and the 'next' step. Future
1072 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1074 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1075 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1076 bCalcVir = bCalcEnerStep ||
1077 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1081 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1082 bCalcVir = bCalcEnerStep ||
1083 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1085 bCalcEner = bCalcEnerStep;
1087 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1089 if (do_ene || do_log || bDoReplEx)
1095 /* Do we need global communication ? */
1096 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1097 do_per_step(step, nstglobalcomm) ||
1098 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1100 force_flags = (GMX_FORCE_STATECHANGED |
1101 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1102 GMX_FORCE_ALLFORCES |
1103 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1104 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1105 (bDoFEP ? GMX_FORCE_DHDL : 0)
1110 /* Now is the time to relax the shells */
1111 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
1112 enforcedRotation, step,
1113 ir, bNS, force_flags, top,
1115 state, f, force_vir, mdatoms,
1116 nrnb, wcycle, graph, groups,
1117 shellfc, fr, t, mu_tot,
1119 ddOpenBalanceRegion, ddCloseBalanceRegion);
1123 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
1124 (or the AWH update will be performed twice for one step when continuing). It would be best to
1125 call this update function from do_md_trajectory_writing but that would occur after do_force.
1126 One would have to divide the update_awh function into one function applying the AWH force
1127 and one doing the AWH bias update. The update AWH bias function could then be called after
1128 do_md_trajectory_writing (then containing update_awh_history).
1129 The checkpointing will in the future probably moved to the start of the md loop which will
1130 rid of this issue. */
1131 if (awh && bCPT && MASTER(cr))
1133 awh->updateHistory(state_global->awhHistory.get());
1136 /* The coordinates (x) are shifted (to get whole molecules)
1138 * This is parallellized as well, and does communication too.
1139 * Check comments in sim_util.c
1141 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation,
1142 step, nrnb, wcycle, top, groups,
1143 state->box, state->x, &state->hist,
1144 f, force_vir, mdatoms, enerd, fcd,
1145 state->lambda, graph,
1146 fr, vsite, mu_tot, t, ed,
1147 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1148 ddOpenBalanceRegion, ddCloseBalanceRegion);
1151 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1152 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1154 rvec *vbuf = nullptr;
1156 wallcycle_start(wcycle, ewcUPDATE);
1157 if (ir->eI == eiVV && bInitStep)
1159 /* if using velocity verlet with full time step Ekin,
1160 * take the first half step only to compute the
1161 * virial for the first step. From there,
1162 * revert back to the initial coordinates
1163 * so that the input is actually the initial step.
1165 snew(vbuf, state->natoms);
1166 copy_rvecn(as_rvec_array(state->v.data()), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1170 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1171 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1174 update_coords(step, ir, mdatoms, state, f, fcd,
1175 ekind, M, upd, etrtVELOCITY1,
1178 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1180 wallcycle_stop(wcycle, ewcUPDATE);
1181 constrain_velocities(step, nullptr,
1185 bCalcVir, do_log, do_ene);
1186 wallcycle_start(wcycle, ewcUPDATE);
1190 /* Need to unshift here if a do_force has been
1191 called in the previous step */
1192 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1194 /* if VV, compute the pressure and constraints */
1195 /* For VV2, we strictly only need this if using pressure
1196 * control, but we really would like to have accurate pressures
1198 * Think about ways around this in the future?
1199 * For now, keep this choice in comments.
1201 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1202 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1204 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1205 if (bCalcEner && ir->eI == eiVVAK)
1207 bSumEkinhOld = TRUE;
1209 /* for vv, the first half of the integration actually corresponds to the previous step.
1210 So we need information from the last step in the first half of the integration */
1211 if (bGStat || do_per_step(step-1, nstglobalcomm))
1213 wallcycle_stop(wcycle, ewcUPDATE);
1214 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1215 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1216 constr, &nullSignaller, state->box,
1217 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1218 (bGStat ? CGLO_GSTAT : 0)
1220 | (bTemp ? CGLO_TEMPERATURE : 0)
1221 | (bPres ? CGLO_PRESSURE : 0)
1222 | (bPres ? CGLO_CONSTRAINT : 0)
1223 | (bStopCM ? CGLO_STOPCM : 0)
1224 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1227 /* explanation of above:
1228 a) We compute Ekin at the full time step
1229 if 1) we are using the AveVel Ekin, and it's not the
1230 initial step, or 2) if we are using AveEkin, but need the full
1231 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1232 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1233 EkinAveVel because it's needed for the pressure */
1234 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1235 top_global, top, state,
1236 &shouldCheckNumberOfBondedInteractions);
1237 wallcycle_start(wcycle, ewcUPDATE);
1239 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1244 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1245 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1247 /* TODO This is only needed when we're about to write
1248 * a checkpoint, because we use it after the restart
1249 * (in a kludge?). But what should we be doing if
1250 * startingFromCheckpoint or bInitStep are true? */
1251 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1253 copy_mat(shake_vir, state->svir_prev);
1254 copy_mat(force_vir, state->fvir_prev);
1256 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1258 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1259 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1260 enerd->term[F_EKIN] = trace(ekind->ekin);
1263 else if (bExchanged)
1265 wallcycle_stop(wcycle, ewcUPDATE);
1266 /* We need the kinetic energy at minus the half step for determining
1267 * the full step kinetic energy and possibly for T-coupling.*/
1268 /* This may not be quite working correctly yet . . . . */
1269 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1270 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1271 constr, &nullSignaller, state->box,
1272 nullptr, &bSumEkinhOld,
1273 CGLO_GSTAT | CGLO_TEMPERATURE);
1274 wallcycle_start(wcycle, ewcUPDATE);
1277 /* if it's the initial step, we performed this first step just to get the constraint virial */
1278 if (ir->eI == eiVV && bInitStep)
1280 copy_rvecn(vbuf, as_rvec_array(state->v.data()), 0, state->natoms);
1283 wallcycle_stop(wcycle, ewcUPDATE);
1286 /* compute the conserved quantity */
1289 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1292 last_ekin = enerd->term[F_EKIN];
1294 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1296 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1298 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1299 if (ir->efep != efepNO && !bRerunMD)
1301 sum_dhdl(enerd, state->lambda, ir->fepvals);
1305 /* ######## END FIRST UPDATE STEP ############## */
1306 /* ######## If doing VV, we now have v(dt) ###### */
1309 /* perform extended ensemble sampling in lambda - we don't
1310 actually move to the new state before outputting
1311 statistics, but if performing simulated tempering, we
1312 do update the velocities and the tau_t. */
1314 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, as_rvec_array(state->v.data()), mdatoms);
1315 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1318 copy_df_history(state_global->dfhist, state->dfhist);
1322 /* Now we have the energies and forces corresponding to the
1323 * coordinates at time t. We must output all of this before
1326 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1327 ir, state, state_global, observablesHistory,
1329 outf, mdebin, ekind, f,
1331 bCPT, bRerunMD, bLastStep,
1332 mdrunOptions.writeConfout,
1334 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1335 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, as_rvec_array(state->x.data()), ir, t, wcycle);
1337 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1338 if (startingFromCheckpoint && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1340 copy_mat(state->svir_prev, shake_vir);
1341 copy_mat(state->fvir_prev, force_vir);
1344 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1346 /* Check whether everything is still allright */
1347 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1353 int nsteps_stop = -1;
1355 /* this just makes signals[].sig compatible with the hack
1356 of sending signals around by MPI_Reduce together with
1358 if ((gmx_get_stop_condition() == gmx_stop_cond_next_ns) ||
1359 (mdrunOptions.reproducible &&
1360 gmx_get_stop_condition() == gmx_stop_cond_next))
1362 /* We need at least two global communication steps to pass
1363 * around the signal. We stop at a pair-list creation step
1364 * to allow for exact continuation, when possible.
1366 signals[eglsSTOPCOND].sig = 1;
1367 nsteps_stop = std::max(ir->nstlist, 2*nstSignalComm);
1369 else if (gmx_get_stop_condition() == gmx_stop_cond_next)
1371 /* Stop directly after the next global communication step.
1372 * This breaks exact continuation.
1374 signals[eglsSTOPCOND].sig = -1;
1375 nsteps_stop = nstSignalComm + 1;
1380 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1381 gmx_get_signal_name(), nsteps_stop);
1385 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1386 gmx_get_signal_name(), nsteps_stop);
1388 handled_stop_condition = (int)gmx_get_stop_condition();
1390 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1391 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1392 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1394 /* Signal to terminate the run */
1395 signals[eglsSTOPCOND].sig = 1;
1398 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1400 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1403 if (bResetCountersHalfMaxH && MASTER(cr) &&
1404 elapsed_time > max_hours*60.0*60.0*0.495)
1406 /* Set flag that will communicate the signal to all ranks in the simulation */
1407 signals[eglsRESETCOUNTERS].sig = 1;
1410 /* In parallel we only have to check for checkpointing in steps
1411 * where we do global communication,
1412 * otherwise the other nodes don't know.
1414 const real cpt_period = mdrunOptions.checkpointOptions.period;
1415 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1418 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1419 signals[eglsCHKPT].set == 0)
1421 signals[eglsCHKPT].sig = 1;
1424 /* ######### START SECOND UPDATE STEP ################# */
1426 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1429 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1431 gmx_bool bIfRandomize;
1432 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1433 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1434 if (constr && bIfRandomize)
1436 constrain_velocities(step, nullptr,
1440 bCalcVir, do_log, do_ene);
1443 /* Box is changed in update() when we do pressure coupling,
1444 * but we should still use the old box for energy corrections and when
1445 * writing it to the energy file, so it matches the trajectory files for
1446 * the same timestep above. Make a copy in a separate array.
1448 copy_mat(state->box, lastbox);
1452 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1454 wallcycle_start(wcycle, ewcUPDATE);
1455 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1458 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1459 /* We can only do Berendsen coupling after we have summed
1460 * the kinetic energy or virial. Since the happens
1461 * in global_state after update, we should only do it at
1462 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1467 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1468 update_pcouple_before_coordinates(fplog, step, ir, state,
1469 parrinellorahmanMu, M,
1475 /* velocity half-step update */
1476 update_coords(step, ir, mdatoms, state, f, fcd,
1477 ekind, M, upd, etrtVELOCITY2,
1481 /* Above, initialize just copies ekinh into ekin,
1482 * it doesn't copy position (for VV),
1483 * and entire integrator for MD.
1486 if (ir->eI == eiVVAK)
1488 /* We probably only need md->homenr, not state->natoms */
1489 if (state->natoms > cbuf_nalloc)
1491 cbuf_nalloc = state->natoms;
1492 srenew(cbuf, cbuf_nalloc);
1494 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1497 update_coords(step, ir, mdatoms, state, f, fcd,
1498 ekind, M, upd, etrtPOSITION, cr, constr);
1499 wallcycle_stop(wcycle, ewcUPDATE);
1501 constrain_coordinates(step, &dvdl_constr, state,
1503 wcycle, upd, constr,
1504 bCalcVir, do_log, do_ene);
1505 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state,
1506 cr, nrnb, wcycle, upd, constr, do_log, do_ene);
1507 finish_update(ir, mdatoms,
1509 nrnb, wcycle, upd, constr);
1511 if (ir->eI == eiVVAK)
1513 /* erase F_EKIN and F_TEMP here? */
1514 /* just compute the kinetic energy at the half step to perform a trotter step */
1515 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1516 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1517 constr, &nullSignaller, lastbox,
1518 nullptr, &bSumEkinhOld,
1519 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1521 wallcycle_start(wcycle, ewcUPDATE);
1522 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1523 /* now we know the scaling, we can compute the positions again again */
1524 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1526 update_coords(step, ir, mdatoms, state, f, fcd,
1527 ekind, M, upd, etrtPOSITION, cr, constr);
1528 wallcycle_stop(wcycle, ewcUPDATE);
1530 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1531 /* are the small terms in the shake_vir here due
1532 * to numerical errors, or are they important
1533 * physically? I'm thinking they are just errors, but not completely sure.
1534 * For now, will call without actually constraining, constr=NULL*/
1535 finish_update(ir, mdatoms,
1537 nrnb, wcycle, upd, nullptr);
1541 /* this factor or 2 correction is necessary
1542 because half of the constraint force is removed
1543 in the vv step, so we have to double it. See
1544 the Redmine issue #1255. It is not yet clear
1545 if the factor of 2 is exact, or just a very
1546 good approximation, and this will be
1547 investigated. The next step is to see if this
1548 can be done adding a dhdl contribution from the
1549 rattle step, but this is somewhat more
1550 complicated with the current code. Will be
1551 investigated, hopefully for 4.6.3. However,
1552 this current solution is much better than
1553 having it completely wrong.
1555 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1559 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1564 /* Need to unshift here */
1565 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1568 if (vsite != nullptr)
1570 wallcycle_start(wcycle, ewcVSITECONSTR);
1571 if (graph != nullptr)
1573 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1575 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
1576 top->idef.iparams, top->idef.il,
1577 fr->ePBC, fr->bMolPBC, cr, state->box);
1579 if (graph != nullptr)
1581 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1583 wallcycle_stop(wcycle, ewcVSITECONSTR);
1586 /* ############## IF NOT VV, Calculate globals HERE ############ */
1587 /* With Leap-Frog we can skip compute_globals at
1588 * non-communication steps, but we need to calculate
1589 * the kinetic energy one step before communication.
1592 // Organize to do inter-simulation signalling on steps if
1593 // and when algorithms require it.
1594 bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1596 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1598 // Since we're already communicating at this step, we
1599 // can propagate intra-simulation signals. Note that
1600 // check_nstglobalcomm has the responsibility for
1601 // choosing the value of nstglobalcomm that is one way
1602 // bGStat becomes true, so we can't get into a
1603 // situation where e.g. checkpointing can't be
1605 bool doIntraSimSignal = true;
1606 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1608 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1609 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1612 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1613 (bGStat ? CGLO_GSTAT : 0)
1614 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1615 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1616 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1617 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1619 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1621 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1622 top_global, top, state,
1623 &shouldCheckNumberOfBondedInteractions);
1627 /* ############# END CALC EKIN AND PRESSURE ################# */
1629 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1630 the virial that should probably be addressed eventually. state->veta has better properies,
1631 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1632 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1634 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1636 /* Sum up the foreign energy and dhdl terms for md and sd.
1637 Currently done every step so that dhdl is correct in the .edr */
1638 sum_dhdl(enerd, state->lambda, ir->fepvals);
1641 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1642 pres, force_vir, shake_vir,
1646 /* ################# END UPDATE STEP 2 ################# */
1647 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1649 /* The coordinates (x) were unshifted in update */
1652 /* We will not sum ekinh_old,
1653 * so signal that we still have to do it.
1655 bSumEkinhOld = TRUE;
1660 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1662 /* use the directly determined last velocity, not actually the averaged half steps */
1663 if (bTrotter && ir->eI == eiVV)
1665 enerd->term[F_EKIN] = last_ekin;
1667 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1669 if (integratorHasConservedEnergyQuantity(ir))
1673 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1677 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1680 /* ######### END PREPARING EDR OUTPUT ########### */
1686 if (fplog && do_log && bDoExpanded)
1688 /* only needed if doing expanded ensemble */
1689 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1690 state_global->dfhist, state->fep_state, ir->nstlog, step);
1694 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1695 t, mdatoms->tmass, enerd, state,
1696 ir->fepvals, ir->expandedvals, lastbox,
1697 shake_vir, force_vir, total_vir, pres,
1698 ekind, mu_tot, constr);
1702 upd_mdebin_step(mdebin);
1705 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1706 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1708 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1710 eprNORMAL, mdebin, fcd, groups, &(ir->opts), awh.get());
1714 pull_print_output(ir->pull_work, step, t);
1717 if (do_per_step(step, ir->nstlog))
1719 if (fflush(fplog) != 0)
1721 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1727 /* Have to do this part _after_ outputting the logfile and the edr file */
1728 /* Gets written into the state at the beginning of next loop*/
1729 state->fep_state = lamnew;
1731 /* Print the remaining wall clock time for the run */
1732 if (isMasterSimMasterRank(ms, cr) &&
1733 (do_verbose || gmx_got_usr_signal()) &&
1738 fprintf(stderr, "\n");
1740 print_time(stderr, walltime_accounting, step, ir, cr);
1743 /* Ion/water position swapping.
1744 * Not done in last step since trajectory writing happens before this call
1745 * in the MD loop and exchanges would be lost anyway. */
1746 bNeedRepartition = FALSE;
1747 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1748 do_per_step(step, ir->swap->nstswap))
1750 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1751 bRerunMD ? rerun_fr.x : as_rvec_array(state->x.data()),
1752 bRerunMD ? rerun_fr.box : state->box,
1753 MASTER(cr) && mdrunOptions.verbose,
1756 if (bNeedRepartition && DOMAINDECOMP(cr))
1758 dd_collect_state(cr->dd, state, state_global);
1762 /* Replica exchange */
1766 bExchanged = replica_exchange(fplog, cr, ms, repl_ex,
1767 state_global, enerd,
1771 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1773 dd_partition_system(fplog, step, cr, TRUE, 1,
1774 state_global, top_global, ir,
1776 state, &f, mdAtoms, top, fr,
1778 nrnb, wcycle, FALSE);
1779 shouldCheckNumberOfBondedInteractions = true;
1780 update_realloc(upd, state->natoms);
1785 startingFromCheckpoint = false;
1787 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1788 /* With all integrators, except VV, we need to retain the pressure
1789 * at the current step for coupling at the next step.
1791 if ((state->flags & (1<<estPRES_PREV)) &&
1793 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1795 /* Store the pressure in t_state for pressure coupling
1796 * at the next MD step.
1798 copy_mat(pres, state->pres_prev);
1801 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1803 if ( (membed != nullptr) && (!bLastStep) )
1805 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1812 /* read next frame from input trajectory */
1813 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1818 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1822 cycles = wallcycle_stop(wcycle, ewcSTEP);
1823 if (DOMAINDECOMP(cr) && wcycle)
1825 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1828 if (!bRerunMD || !rerun_fr.bStep)
1830 /* increase the MD step number */
1835 /* TODO make a counter-reset module */
1836 /* If it is time to reset counters, set a flag that remains
1837 true until counters actually get reset */
1838 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1839 signals[eglsRESETCOUNTERS].set != 0)
1841 if (pme_loadbal_is_active(pme_loadbal))
1843 /* Do not permit counter reset while PME load
1844 * balancing is active. The only purpose for resetting
1845 * counters is to measure reliable performance data,
1846 * and that can't be done before balancing
1849 * TODO consider fixing this by delaying the reset
1850 * until after load balancing completes,
1851 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1852 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1853 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1854 "resetting counters later in the run, e.g. with gmx "
1855 "mdrun -resetstep.", step);
1857 reset_all_counters(fplog, mdlog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1858 use_GPU(fr->nbv) ? fr->nbv : nullptr, fr->pmedata);
1859 wcycle_set_reset_counters(wcycle, -1);
1860 if (!thisRankHasDuty(cr, DUTY_PME))
1862 /* Tell our PME node to reset its counters */
1863 gmx_pme_send_resetcounters(cr, step);
1865 /* Correct max_hours for the elapsed time */
1866 max_hours -= elapsed_time/(60.0*60.0);
1867 /* If mdrun -maxh -resethway was active, it can only trigger once */
1868 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1869 /* Reset can only happen once, so clear the triggering flag. */
1870 signals[eglsRESETCOUNTERS].set = 0;
1873 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1874 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1877 /* End of main MD loop */
1879 /* Closing TNG files can include compressing data. Therefore it is good to do that
1880 * before stopping the time measurements. */
1881 mdoutf_tng_close(outf);
1883 /* Stop measuring walltime */
1884 walltime_accounting_end(walltime_accounting);
1886 if (bRerunMD && MASTER(cr))
1891 if (!thisRankHasDuty(cr, DUTY_PME))
1893 /* Tell the PME only node to finish */
1894 gmx_pme_send_finish(cr);
1899 if (ir->nstcalcenergy > 0 && !bRerunMD)
1901 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1902 eprAVER, mdebin, fcd, groups, &(ir->opts), awh.get());
1905 done_mdebin(mdebin);
1910 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
1913 done_shellfc(fplog, shellfc, step_rel);
1915 if (useReplicaExchange && MASTER(cr))
1917 print_replica_exchange_statistics(fplog, repl_ex);
1920 // Clean up swapcoords
1921 if (ir->eSwapCoords != eswapNO)
1923 finish_swapcoords(ir->swap);
1926 /* Do essential dynamics cleanup if needed. Close .edo file */
1929 /* IMD cleanup, if bIMD is TRUE. */
1930 IMD_finalize(ir->bIMD, ir->imd);
1932 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1933 if (step_rel >= wcycle_get_reset_counters(wcycle) &&
1934 signals[eglsRESETCOUNTERS].set == 0 &&
1935 !bResetCountersHalfMaxH)
1937 walltime_accounting_set_valid_finish(walltime_accounting);
1940 destroy_enerdata(enerd);