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
55 #include "gromacs/awh/awh.h"
56 #include "gromacs/commandline/filenm.h"
57 #include "gromacs/compat/make_unique.h"
58 #include "gromacs/domdec/collect.h"
59 #include "gromacs/domdec/domdec.h"
60 #include "gromacs/domdec/domdec_network.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/essentialdynamics/edsam.h"
63 #include "gromacs/ewald/pme.h"
64 #include "gromacs/ewald/pme-load-balancing.h"
65 #include "gromacs/fileio/trxio.h"
66 #include "gromacs/gmxlib/network.h"
67 #include "gromacs/gmxlib/nrnb.h"
68 #include "gromacs/gpu_utils/gpu_utils.h"
69 #include "gromacs/imd/imd.h"
70 #include "gromacs/listed-forces/manage-threading.h"
71 #include "gromacs/math/functions.h"
72 #include "gromacs/math/utilities.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/math/vectypes.h"
75 #include "gromacs/mdlib/compute_io.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/ebin.h"
78 #include "gromacs/mdlib/expanded.h"
79 #include "gromacs/mdlib/force.h"
80 #include "gromacs/mdlib/force_flags.h"
81 #include "gromacs/mdlib/forcerec.h"
82 #include "gromacs/mdlib/md_support.h"
83 #include "gromacs/mdlib/mdatoms.h"
84 #include "gromacs/mdlib/mdebin.h"
85 #include "gromacs/mdlib/mdoutf.h"
86 #include "gromacs/mdlib/mdrun.h"
87 #include "gromacs/mdlib/mdsetup.h"
88 #include "gromacs/mdlib/membed.h"
89 #include "gromacs/mdlib/nb_verlet.h"
90 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
91 #include "gromacs/mdlib/ns.h"
92 #include "gromacs/mdlib/repl_ex.h"
93 #include "gromacs/mdlib/shellfc.h"
94 #include "gromacs/mdlib/sighandler.h"
95 #include "gromacs/mdlib/sim_util.h"
96 #include "gromacs/mdlib/simulationsignal.h"
97 #include "gromacs/mdlib/tgroup.h"
98 #include "gromacs/mdlib/trajectory_writing.h"
99 #include "gromacs/mdlib/update.h"
100 #include "gromacs/mdlib/vcm.h"
101 #include "gromacs/mdlib/vsite.h"
102 #include "gromacs/mdtypes/awh-history.h"
103 #include "gromacs/mdtypes/awh-params.h"
104 #include "gromacs/mdtypes/commrec.h"
105 #include "gromacs/mdtypes/df_history.h"
106 #include "gromacs/mdtypes/energyhistory.h"
107 #include "gromacs/mdtypes/fcdata.h"
108 #include "gromacs/mdtypes/forcerec.h"
109 #include "gromacs/mdtypes/group.h"
110 #include "gromacs/mdtypes/inputrec.h"
111 #include "gromacs/mdtypes/interaction_const.h"
112 #include "gromacs/mdtypes/md_enums.h"
113 #include "gromacs/mdtypes/mdatom.h"
114 #include "gromacs/mdtypes/observableshistory.h"
115 #include "gromacs/mdtypes/state.h"
116 #include "gromacs/pbcutil/mshift.h"
117 #include "gromacs/pbcutil/pbc.h"
118 #include "gromacs/pulling/pull.h"
119 #include "gromacs/swap/swapcoords.h"
120 #include "gromacs/timing/wallcycle.h"
121 #include "gromacs/timing/walltime_accounting.h"
122 #include "gromacs/topology/atoms.h"
123 #include "gromacs/topology/idef.h"
124 #include "gromacs/topology/mtop_util.h"
125 #include "gromacs/topology/topology.h"
126 #include "gromacs/trajectory/trajectoryframe.h"
127 #include "gromacs/utility/basedefinitions.h"
128 #include "gromacs/utility/cstringutil.h"
129 #include "gromacs/utility/fatalerror.h"
130 #include "gromacs/utility/logger.h"
131 #include "gromacs/utility/real.h"
132 #include "gromacs/utility/smalloc.h"
134 #include "integrator.h"
137 #include "corewrap.h"
140 using gmx::SimulationSignaller;
142 //! Resets all the counters.
143 static void reset_all_counters(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
145 int64_t *step_rel, t_inputrec *ir,
146 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
147 gmx_walltime_accounting_t walltime_accounting,
148 struct nonbonded_verlet_t *nbv,
149 struct gmx_pme_t *pme)
151 char sbuf[STEPSTRSIZE];
153 /* Reset all the counters related to performance over the run */
154 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
155 "step %s: resetting all time and cycle counters",
156 gmx_step_str(step, sbuf));
160 nbnxn_gpu_reset_timings(nbv);
163 if (pme_gpu_task_enabled(pme))
165 pme_gpu_reset_timings(pme);
168 if (use_GPU(nbv) || pme_gpu_task_enabled(pme))
173 wallcycle_stop(wcycle, ewcRUN);
174 wallcycle_reset_all(wcycle);
175 if (DOMAINDECOMP(cr))
177 reset_dd_statistics_counters(cr->dd);
180 ir->init_step += *step_rel;
181 ir->nsteps -= *step_rel;
183 wallcycle_start(wcycle, ewcRUN);
184 walltime_accounting_start(walltime_accounting);
185 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
188 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
190 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
191 * \param[in,out] globalState The global state container
192 * \param[in] constructVsites When true, vsite coordinates are constructed
193 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
194 * \param[in] idef Topology parameters, used for constructing vsites
195 * \param[in] timeStep Time step, used for constructing vsites
196 * \param[in] forceRec Force record, used for constructing vsites
197 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
198 * \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
200 static void prepareRerunState(const t_trxframe &rerunFrame,
201 t_state *globalState,
202 bool constructVsites,
203 const gmx_vsite_t *vsite,
206 const t_forcerec &forceRec,
208 gmx_bool *warnWhenNoV)
210 for (int i = 0; i < globalState->natoms; i++)
212 copy_rvec(rerunFrame.x[i], globalState->x[i]);
216 for (int i = 0; i < globalState->natoms; i++)
218 copy_rvec(rerunFrame.v[i], globalState->v[i]);
223 for (int i = 0; i < globalState->natoms; i++)
225 clear_rvec(globalState->v[i]);
229 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
230 " Ekin, temperature and pressure are incorrect,\n"
231 " the virial will be incorrect when constraints are present.\n"
233 *warnWhenNoV = FALSE;
236 copy_mat(rerunFrame.box, globalState->box);
240 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
244 /* Following is necessary because the graph may get out of sync
245 * with the coordinates if we only have every N'th coordinate set
247 mk_mshift(nullptr, graph, forceRec.ePBC, globalState->box, as_rvec_array(globalState->x.data()));
248 shift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
250 construct_vsites(vsite, as_rvec_array(globalState->x.data()), timeStep, as_rvec_array(globalState->v.data()),
251 idef.iparams, idef.il,
252 forceRec.ePBC, forceRec.bMolPBC, nullptr, globalState->box);
255 unshift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
260 void gmx::Integrator::do_md()
262 // TODO Historically, the EM and MD "integrators" used different
263 // names for the t_inputrec *parameter, but these must have the
264 // same name, now that it's a member of a struct. We use this ir
265 // alias to avoid a large ripple of nearly useless changes.
266 // t_inputrec is being replaced by IMdpOptionsProvider, so this
267 // will go away eventually.
268 t_inputrec *ir = inputrec;
269 gmx_mdoutf *outf = nullptr;
270 int64_t step, step_rel;
272 double t, t0, lam0[efptNR];
273 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
274 gmx_bool bNS, bNStList, bSimAnn, bStopCM,
275 bFirstStep, bInitStep, bLastStep = FALSE;
276 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
277 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
278 bForceUpdate = FALSE, bCPT;
279 gmx_bool bMasterState;
280 int force_flags, cglo_flags;
281 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
286 matrix parrinellorahmanMu, M;
288 gmx_repl_ex_t repl_ex = nullptr;
291 t_mdebin *mdebin = nullptr;
292 gmx_enerdata_t *enerd;
293 PaddedRVecVector f {};
294 gmx_global_stat_t gstat;
295 gmx_update_t *upd = nullptr;
296 t_graph *graph = nullptr;
297 gmx_groups_t *groups;
298 gmx_ekindata_t *ekind;
299 gmx_shellfc_t *shellfc;
300 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
301 gmx_bool bResetCountersHalfMaxH = FALSE;
302 gmx_bool bTemp, bPres, bTrotter;
304 rvec *cbuf = nullptr;
311 real saved_conserved_quantity = 0;
315 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
316 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
319 /* PME load balancing data for GPU kernels */
320 pme_load_balancing_t *pme_loadbal = nullptr;
321 gmx_bool bPMETune = FALSE;
322 gmx_bool bPMETunePrinting = FALSE;
325 gmx_bool bIMDstep = FALSE;
328 /* Temporary addition for FAHCORE checkpointing */
331 /* Domain decomposition could incorrectly miss a bonded
332 interaction, but checking for that requires a global
333 communication stage, which does not otherwise happen in DD
334 code. So we do that alongside the first global energy reduction
335 after a new DD is made. These variables handle whether the
336 check happens, and the result it returns. */
337 bool shouldCheckNumberOfBondedInteractions = false;
338 int totalNumberOfBondedInteractions = -1;
340 SimulationSignals signals;
341 // Most global communnication stages don't propagate mdrun
342 // signals, and will use this object to achieve that.
343 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
345 if (!mdrunOptions.writeConfout)
347 // This is on by default, and the main known use case for
348 // turning it off is for convenience in benchmarking, which is
349 // something that should not show up in the general user
351 GMX_LOG(mdlog.info).asParagraph().
352 appendText("The -noconfout functionality is deprecated, and may be removed in a future version.");
354 if (mdrunOptions.timingOptions.resetHalfway)
356 GMX_LOG(mdlog.info).asParagraph().
357 appendText("The -resethway functionality is deprecated, and may be removed in a future version.");
360 /* Signal to reset the counters half the simulation steps. */
361 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
363 /* Signal to reset the counters halfway the simulation time. */
364 bResetCountersHalfMaxH = (mdrunOptions.maximumHoursToRun > 0);
367 /* md-vv uses averaged full step velocities for T-control
368 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
369 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
370 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
372 const gmx_bool bRerunMD = mdrunOptions.rerun;
373 int nstglobalcomm = mdrunOptions.globalCommunicationInterval;
377 /* Since we don't know if the frames read are related in any way,
378 * rebuild the neighborlist at every step.
381 ir->nstcalcenergy = 1;
385 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
386 bGStatEveryStep = (nstglobalcomm == 1);
390 ir->nstxout_compressed = 0;
392 groups = &top_global->groups;
394 std::unique_ptr<EssentialDynamics> ed = nullptr;
395 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
397 /* Initialize essential dynamics sampling */
398 ed = init_edsam(opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm),
401 state_global, observablesHistory,
402 oenv, mdrunOptions.continuationOptions.appendFiles);
406 init_md(fplog, cr, outputProvider, ir, oenv, mdrunOptions,
407 &t, &t0, state_global, lam0,
408 nrnb, top_global, &upd, deform,
409 nfile, fnm, &outf, &mdebin,
410 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, wcycle);
412 clear_mat(total_vir);
414 /* Energy terms and groups */
416 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
419 /* Kinetic energy data */
421 init_ekindata(fplog, top_global, &(ir->opts), ekind);
422 /* Copy the cos acceleration to the groups struct */
423 ekind->cosacc.cos_accel = ir->cos_accel;
425 gstat = global_stat_init(ir);
427 /* Check for polarizable models and flexible constraints */
428 shellfc = init_shell_flexcon(fplog,
429 top_global, constr ? constr->numFlexibleConstraints() : 0,
430 ir->nstcalcenergy, DOMAINDECOMP(cr));
433 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
434 if ((io > 2000) && MASTER(cr))
437 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
442 /* Set up interactive MD (IMD) */
443 init_IMD(ir, cr, ms, top_global, fplog, ir->nstcalcenergy,
444 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
445 nfile, fnm, oenv, mdrunOptions);
447 // Local state only becomes valid now.
448 std::unique_ptr<t_state> stateInstance;
451 if (DOMAINDECOMP(cr))
453 top = dd_init_local_top(top_global);
455 stateInstance = compat::make_unique<t_state>();
456 state = stateInstance.get();
457 dd_init_local_state(cr->dd, state_global, state);
459 /* Distribute the charge groups over the nodes from the master node */
460 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
461 state_global, top_global, ir,
462 state, &f, mdAtoms, top, fr,
464 nrnb, nullptr, FALSE);
465 shouldCheckNumberOfBondedInteractions = true;
466 update_realloc(upd, state->natoms);
470 state_change_natoms(state_global, state_global->natoms);
471 /* We need to allocate one element extra, since we might use
472 * (unaligned) 4-wide SIMD loads to access rvec entries.
474 f.resize(gmx::paddedRVecVectorSize(state_global->natoms));
475 /* Copy the pointer to the global state */
476 state = state_global;
479 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
480 &graph, mdAtoms, constr, vsite, shellfc);
482 update_realloc(upd, state->natoms);
485 auto mdatoms = mdAtoms->mdatoms();
487 // NOTE: The global state is no longer used at this point.
488 // But state_global is still used as temporary storage space for writing
489 // the global state to file and potentially for replica exchange.
490 // (Global topology should persist.)
492 update_mdatoms(mdatoms, state->lambda[efptMASS]);
494 const ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
495 bool startingFromCheckpoint = continuationOptions.startedFromCheckpoint;
499 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
504 if (startingFromCheckpoint)
506 /* Update mdebin with energy history if appending to output files */
507 if (continuationOptions.appendFiles)
509 restore_energyhistory_from_state(mdebin, observablesHistory->energyHistory.get());
511 else if (observablesHistory->energyHistory != nullptr)
513 /* We might have read an energy history from checkpoint.
514 * As we are not appending, we want to restart the statistics.
515 * Free the allocated memory and reset the counts.
517 observablesHistory->energyHistory = {};
520 if (observablesHistory->energyHistory == nullptr)
522 observablesHistory->energyHistory = compat::make_unique<energyhistory_t>();
524 /* Set the initial energy history in state by updating once */
525 update_energyhistory(observablesHistory->energyHistory.get(), mdebin);
528 // TODO: Remove this by converting AWH into a ForceProvider
529 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms, startingFromCheckpoint,
531 opt2fn("-awh", nfile, fnm), ir->pull_work);
533 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
534 if (useReplicaExchange && MASTER(cr))
536 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir,
539 /* PME tuning is only supported in the Verlet scheme, with PME for
540 * Coulomb. It is not supported with only LJ PME, or for
542 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !bRerunMD &&
543 !mdrunOptions.reproducible && ir->cutoff_scheme != ecutsGROUP);
546 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
547 fr->ic, fr->nbv->listParams.get(), fr->pmedata, use_GPU(fr->nbv),
551 if (!ir->bContinuation && !bRerunMD)
553 if (state->flags & (1 << estV))
555 /* Set the velocities of vsites, shells and frozen atoms to zero */
556 for (i = 0; i < mdatoms->homenr; i++)
558 if (mdatoms->ptype[i] == eptVSite ||
559 mdatoms->ptype[i] == eptShell)
561 clear_rvec(state->v[i]);
563 else if (mdatoms->cFREEZE)
565 for (m = 0; m < DIM; m++)
567 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
578 /* Constrain the initial coordinates and velocities */
579 do_constrain_first(fplog, constr, ir, mdatoms, state);
583 /* Construct the virtual sites for the initial configuration */
584 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, nullptr,
585 top->idef.iparams, top->idef.il,
586 fr->ePBC, fr->bMolPBC, cr, state->box);
590 if (ir->efep != efepNO)
592 /* Set free energy calculation frequency as the greatest common
593 * denominator of nstdhdl and repl_ex_nst. */
594 nstfep = ir->fepvals->nstdhdl;
597 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
599 if (useReplicaExchange)
601 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
605 /* Be REALLY careful about what flags you set here. You CANNOT assume
606 * this is the first step, since we might be restarting from a checkpoint,
607 * and in that case we should not do any modifications to the state.
609 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
611 if (continuationOptions.haveReadEkin)
613 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
616 cglo_flags = (CGLO_INITIALIZATION | CGLO_TEMPERATURE | CGLO_GSTAT
617 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
618 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
619 | (continuationOptions.haveReadEkin ? CGLO_READEKIN : 0));
621 bSumEkinhOld = FALSE;
622 /* To minimize communication, compute_globals computes the COM velocity
623 * and the kinetic energy for the velocities without COM motion removed.
624 * Thus to get the kinetic energy without the COM contribution, we need
625 * to call compute_globals twice.
627 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
629 int cglo_flags_iteration = cglo_flags;
630 if (bStopCM && cgloIteration == 0)
632 cglo_flags_iteration |= CGLO_STOPCM;
633 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
635 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
636 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
637 constr, &nullSignaller, state->box,
638 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags_iteration
639 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
641 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
642 top_global, top, state,
643 &shouldCheckNumberOfBondedInteractions);
644 if (ir->eI == eiVVAK)
646 /* a second call to get the half step temperature initialized as well */
647 /* we do the same call as above, but turn the pressure off -- internally to
648 compute_globals, this is recognized as a velocity verlet half-step
649 kinetic energy calculation. This minimized excess variables, but
650 perhaps loses some logic?*/
652 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
653 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
654 constr, &nullSignaller, state->box,
655 nullptr, &bSumEkinhOld,
656 cglo_flags & ~CGLO_PRESSURE);
659 /* Calculate the initial half step temperature, and save the ekinh_old */
660 if (!continuationOptions.startedFromCheckpoint)
662 for (i = 0; (i < ir->opts.ngtc); i++)
664 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
668 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
669 temperature control */
670 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
674 if (!ir->bContinuation)
676 if (constr && ir->eConstrAlg == econtLINCS)
679 "RMS relative constraint deviation after constraining: %.2e\n",
682 if (EI_STATE_VELOCITY(ir->eI))
684 real temp = enerd->term[F_TEMP];
687 /* Result of Ekin averaged over velocities of -half
688 * and +half step, while we only have -half step here.
692 fprintf(fplog, "Initial temperature: %g K\n", temp);
698 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
699 " input trajectory '%s'\n\n",
700 *(top_global->name), opt2fn("-rerun", nfile, fnm));
701 if (mdrunOptions.verbose)
703 fprintf(stderr, "Calculated time to finish depends on nsteps from "
704 "run input file,\nwhich may not correspond to the time "
705 "needed to process input trajectory.\n\n");
711 fprintf(stderr, "starting mdrun '%s'\n",
712 *(top_global->name));
715 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
719 sprintf(tbuf, "%s", "infinite");
721 if (ir->init_step > 0)
723 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
724 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
725 gmx_step_str(ir->init_step, sbuf2),
726 ir->init_step*ir->delta_t);
730 fprintf(stderr, "%s steps, %s ps.\n",
731 gmx_step_str(ir->nsteps, sbuf), tbuf);
734 fprintf(fplog, "\n");
737 walltime_accounting_start(walltime_accounting);
738 wallcycle_start(wcycle, ewcRUN);
739 print_start(fplog, cr, walltime_accounting, "mdrun");
741 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
743 chkpt_ret = fcCheckPointParallel( cr->nodeid,
747 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
751 /***********************************************************
755 ************************************************************/
757 /* if rerunMD then read coordinates and velocities from input trajectory */
760 if (getenv("GMX_FORCE_UPDATE"))
768 bLastStep = !read_first_frame(oenv, &status,
769 opt2fn("-rerun", nfile, fnm),
770 &rerun_fr, TRX_NEED_X | TRX_READ_V);
771 if (rerun_fr.natoms != top_global->natoms)
774 "Number of atoms in trajectory (%d) does not match the "
775 "run input file (%d)\n",
776 rerun_fr.natoms, top_global->natoms);
778 if (ir->ePBC != epbcNONE)
782 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);
784 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
786 gmx_fatal(FARGS, "Rerun trajectory frame step %ld time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
793 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
796 if (ir->ePBC != epbcNONE)
798 /* Set the shift vectors.
799 * Necessary here when have a static box different from the tpr box.
801 calc_shifts(rerun_fr.box, fr->shift_vec);
805 /* Loop over MD steps or if rerunMD to end of input trajectory,
806 * or, if max_hours>0, until max_hours is reached.
808 real max_hours = mdrunOptions.maximumHoursToRun;
810 /* Skip the first Nose-Hoover integration when we get the state from tpx */
811 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
812 bSumEkinhOld = FALSE;
814 bNeedRepartition = FALSE;
816 bool simulationsShareState = false;
817 int nstSignalComm = nstglobalcomm;
819 // TODO This implementation of ensemble orientation restraints is nasty because
820 // a user can't just do multi-sim with single-sim orientation restraints.
821 bool usingEnsembleRestraints = (fcd->disres.nsystems > 1) || (ms && fcd->orires.nr);
822 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && ms);
824 // Replica exchange, ensemble restraints and AWH need all
825 // simulations to remain synchronized, so they need
826 // checkpoints and stop conditions to act on the same step, so
827 // the propagation of such signals must take place between
828 // simulations, not just within simulations.
829 // TODO: Make algorithm initializers set these flags.
830 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
831 bool resetCountersIsLocal = true;
832 signals[eglsCHKPT] = SimulationSignal(!simulationsShareState);
833 signals[eglsSTOPCOND] = SimulationSignal(!simulationsShareState);
834 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
836 if (simulationsShareState)
838 // Inter-simulation signal communication does not need to happen
839 // often, so we use a minimum of 200 steps to reduce overhead.
840 const int c_minimumInterSimulationSignallingInterval = 200;
841 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1)/nstglobalcomm)*nstglobalcomm;
845 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion = (DOMAINDECOMP(cr) ? DdOpenBalanceRegionBeforeForceComputation::yes : DdOpenBalanceRegionBeforeForceComputation::no);
846 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion = (DOMAINDECOMP(cr) ? DdCloseBalanceRegionAfterForceComputation::yes : DdCloseBalanceRegionAfterForceComputation::no);
848 step = ir->init_step;
851 // TODO extract this to new multi-simulation module
852 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
854 if (!multisim_int_all_are_equal(ms, ir->nsteps))
856 GMX_LOG(mdlog.warning).appendText(
857 "Note: The number of steps is not consistent across multi simulations,\n"
858 "but we are proceeding anyway!");
860 if (!multisim_int_all_are_equal(ms, ir->init_step))
862 GMX_LOG(mdlog.warning).appendText(
863 "Note: The initial step is not consistent across multi simulations,\n"
864 "but we are proceeding anyway!");
868 /* and stop now if we should */
869 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
873 /* Determine if this is a neighbor search step */
874 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
876 if (bPMETune && bNStList)
878 /* PME grid + cut-off optimization with GPUs or PME nodes */
879 pme_loadbal_do(pme_loadbal, cr,
880 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
888 wallcycle_start(wcycle, ewcSTEP);
894 step = rerun_fr.step;
895 step_rel = step - ir->init_step;
908 bLastStep = (step_rel == ir->nsteps);
909 t = t0 + step*ir->delta_t;
912 // TODO Refactor this, so that nstfep does not need a default value of zero
913 if (ir->efep != efepNO || ir->bSimTemp)
915 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
916 requiring different logic. */
921 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
926 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
928 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
929 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
930 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
931 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
934 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep &&
935 do_per_step(step, replExParams.exchangeInterval));
939 update_annealing_target_temp(ir, t, upd);
942 if (bRerunMD && MASTER(cr))
944 const bool constructVsites = (vsite && mdrunOptions.rerunConstructVsites);
945 if (constructVsites && DOMAINDECOMP(cr))
947 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
949 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top->idef, ir->delta_t, *fr, graph, &bRerunWarnNoV);
952 /* Stop Center of Mass motion */
953 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
957 /* for rerun MD always do Neighbour Searching */
958 bNS = (bFirstStep || ir->nstlist != 0);
962 /* Determine whether or not to do Neighbour Searching */
963 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
966 /* < 0 means stop at next step, > 0 means stop at next NS step */
967 if ( (signals[eglsSTOPCOND].set < 0) ||
968 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
973 /* do_log triggers energy and virial calculation. Because this leads
974 * to different code paths, forces can be different. Thus for exact
975 * continuation we should avoid extra log output.
976 * Note that the || bLastStep can result in non-exact continuation
977 * beyond the last step. But we don't consider that to be an issue.
979 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
980 do_verbose = mdrunOptions.verbose &&
981 (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep || bRerunMD);
983 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
991 bMasterState = FALSE;
992 /* Correct the new box if it is too skewed */
993 if (inputrecDynamicBox(ir))
995 if (correct_box(fplog, step, state->box, graph))
1000 if (DOMAINDECOMP(cr) && bMasterState)
1002 dd_collect_state(cr->dd, state, state_global);
1006 if (DOMAINDECOMP(cr))
1008 /* Repartition the domain decomposition */
1009 dd_partition_system(fplog, step, cr,
1010 bMasterState, nstglobalcomm,
1011 state_global, top_global, ir,
1012 state, &f, mdAtoms, top, fr,
1015 do_verbose && !bPMETunePrinting);
1016 shouldCheckNumberOfBondedInteractions = true;
1017 update_realloc(upd, state->natoms);
1021 if (MASTER(cr) && do_log)
1023 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
1026 if (ir->efep != efepNO)
1028 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1031 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1034 /* We need the kinetic energy at minus the half step for determining
1035 * the full step kinetic energy and possibly for T-coupling.*/
1036 /* This may not be quite working correctly yet . . . . */
1037 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1038 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1039 constr, &nullSignaller, state->box,
1040 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1041 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1042 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1043 top_global, top, state,
1044 &shouldCheckNumberOfBondedInteractions);
1046 clear_mat(force_vir);
1048 /* We write a checkpoint at this MD step when:
1049 * either at an NS step when we signalled through gs,
1050 * or at the last step (but not when we do not want confout),
1051 * but never at the first step or with rerun.
1053 bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1054 (bLastStep && mdrunOptions.writeConfout)) &&
1055 step > ir->init_step && !bRerunMD);
1058 signals[eglsCHKPT].set = 0;
1061 /* Determine the energy and pressure:
1062 * at nstcalcenergy steps and at energy output steps (set below).
1064 if (EI_VV(ir->eI) && (!bInitStep))
1066 /* for vv, the first half of the integration actually corresponds
1067 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1068 but the virial needs to be calculated on both the current step and the 'next' step. Future
1069 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1071 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1072 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1073 bCalcVir = bCalcEnerStep ||
1074 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1078 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1079 bCalcVir = bCalcEnerStep ||
1080 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1082 bCalcEner = bCalcEnerStep;
1084 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1086 if (do_ene || do_log || bDoReplEx)
1092 /* Do we need global communication ? */
1093 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1094 do_per_step(step, nstglobalcomm) ||
1095 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1097 force_flags = (GMX_FORCE_STATECHANGED |
1098 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1099 GMX_FORCE_ALLFORCES |
1100 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1101 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1102 (bDoFEP ? GMX_FORCE_DHDL : 0)
1107 /* Now is the time to relax the shells */
1108 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
1109 enforcedRotation, step,
1110 ir, bNS, force_flags, top,
1112 state, f, force_vir, mdatoms,
1113 nrnb, wcycle, graph, groups,
1114 shellfc, fr, t, mu_tot,
1116 ddOpenBalanceRegion, ddCloseBalanceRegion);
1120 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
1121 (or the AWH update will be performed twice for one step when continuing). It would be best to
1122 call this update function from do_md_trajectory_writing but that would occur after do_force.
1123 One would have to divide the update_awh function into one function applying the AWH force
1124 and one doing the AWH bias update. The update AWH bias function could then be called after
1125 do_md_trajectory_writing (then containing update_awh_history).
1126 The checkpointing will in the future probably moved to the start of the md loop which will
1127 rid of this issue. */
1128 if (awh && bCPT && MASTER(cr))
1130 awh->updateHistory(state_global->awhHistory.get());
1133 /* The coordinates (x) are shifted (to get whole molecules)
1135 * This is parallellized as well, and does communication too.
1136 * Check comments in sim_util.c
1138 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation,
1139 step, nrnb, wcycle, top, groups,
1140 state->box, state->x, &state->hist,
1141 f, force_vir, mdatoms, enerd, fcd,
1142 state->lambda, graph,
1143 fr, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
1144 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1145 ddOpenBalanceRegion, ddCloseBalanceRegion);
1148 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1149 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1151 rvec *vbuf = nullptr;
1153 wallcycle_start(wcycle, ewcUPDATE);
1154 if (ir->eI == eiVV && bInitStep)
1156 /* if using velocity verlet with full time step Ekin,
1157 * take the first half step only to compute the
1158 * virial for the first step. From there,
1159 * revert back to the initial coordinates
1160 * so that the input is actually the initial step.
1162 snew(vbuf, state->natoms);
1163 copy_rvecn(as_rvec_array(state->v.data()), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1167 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1168 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1171 update_coords(step, ir, mdatoms, state, f, fcd,
1172 ekind, M, upd, etrtVELOCITY1,
1175 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1177 wallcycle_stop(wcycle, ewcUPDATE);
1178 constrain_velocities(step, nullptr,
1182 bCalcVir, do_log, do_ene);
1183 wallcycle_start(wcycle, ewcUPDATE);
1187 /* Need to unshift here if a do_force has been
1188 called in the previous step */
1189 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1191 /* if VV, compute the pressure and constraints */
1192 /* For VV2, we strictly only need this if using pressure
1193 * control, but we really would like to have accurate pressures
1195 * Think about ways around this in the future?
1196 * For now, keep this choice in comments.
1198 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1199 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1201 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1202 if (bCalcEner && ir->eI == eiVVAK)
1204 bSumEkinhOld = TRUE;
1206 /* for vv, the first half of the integration actually corresponds to the previous step.
1207 So we need information from the last step in the first half of the integration */
1208 if (bGStat || do_per_step(step-1, nstglobalcomm))
1210 wallcycle_stop(wcycle, ewcUPDATE);
1211 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1212 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1213 constr, &nullSignaller, state->box,
1214 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1215 (bGStat ? CGLO_GSTAT : 0)
1217 | (bTemp ? CGLO_TEMPERATURE : 0)
1218 | (bPres ? CGLO_PRESSURE : 0)
1219 | (bPres ? CGLO_CONSTRAINT : 0)
1220 | (bStopCM ? CGLO_STOPCM : 0)
1221 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1224 /* explanation of above:
1225 a) We compute Ekin at the full time step
1226 if 1) we are using the AveVel Ekin, and it's not the
1227 initial step, or 2) if we are using AveEkin, but need the full
1228 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1229 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1230 EkinAveVel because it's needed for the pressure */
1231 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1232 top_global, top, state,
1233 &shouldCheckNumberOfBondedInteractions);
1234 wallcycle_start(wcycle, ewcUPDATE);
1236 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1241 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1242 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1244 /* TODO This is only needed when we're about to write
1245 * a checkpoint, because we use it after the restart
1246 * (in a kludge?). But what should we be doing if
1247 * startingFromCheckpoint or bInitStep are true? */
1248 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1250 copy_mat(shake_vir, state->svir_prev);
1251 copy_mat(force_vir, state->fvir_prev);
1253 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1255 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1256 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1257 enerd->term[F_EKIN] = trace(ekind->ekin);
1260 else if (bExchanged)
1262 wallcycle_stop(wcycle, ewcUPDATE);
1263 /* We need the kinetic energy at minus the half step for determining
1264 * the full step kinetic energy and possibly for T-coupling.*/
1265 /* This may not be quite working correctly yet . . . . */
1266 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1267 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1268 constr, &nullSignaller, state->box,
1269 nullptr, &bSumEkinhOld,
1270 CGLO_GSTAT | CGLO_TEMPERATURE);
1271 wallcycle_start(wcycle, ewcUPDATE);
1274 /* if it's the initial step, we performed this first step just to get the constraint virial */
1275 if (ir->eI == eiVV && bInitStep)
1277 copy_rvecn(vbuf, as_rvec_array(state->v.data()), 0, state->natoms);
1280 wallcycle_stop(wcycle, ewcUPDATE);
1283 /* compute the conserved quantity */
1286 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1289 last_ekin = enerd->term[F_EKIN];
1291 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1293 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1295 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1296 if (ir->efep != efepNO && !bRerunMD)
1298 sum_dhdl(enerd, state->lambda, ir->fepvals);
1302 /* ######## END FIRST UPDATE STEP ############## */
1303 /* ######## If doing VV, we now have v(dt) ###### */
1306 /* perform extended ensemble sampling in lambda - we don't
1307 actually move to the new state before outputting
1308 statistics, but if performing simulated tempering, we
1309 do update the velocities and the tau_t. */
1311 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, as_rvec_array(state->v.data()), mdatoms);
1312 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1315 copy_df_history(state_global->dfhist, state->dfhist);
1319 /* Now we have the energies and forces corresponding to the
1320 * coordinates at time t. We must output all of this before
1323 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1324 ir, state, state_global, observablesHistory,
1326 outf, mdebin, ekind, f,
1328 bCPT, bRerunMD, bLastStep,
1329 mdrunOptions.writeConfout,
1331 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1332 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, as_rvec_array(state->x.data()), ir, t, wcycle);
1334 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1335 if (startingFromCheckpoint && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1337 copy_mat(state->svir_prev, shake_vir);
1338 copy_mat(state->fvir_prev, force_vir);
1341 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1343 /* Check whether everything is still allright */
1344 if ((static_cast<int>(gmx_get_stop_condition()) > handled_stop_condition)
1350 int nsteps_stop = -1;
1352 /* this just makes signals[].sig compatible with the hack
1353 of sending signals around by MPI_Reduce together with
1355 if ((gmx_get_stop_condition() == gmx_stop_cond_next_ns) ||
1356 (mdrunOptions.reproducible &&
1357 gmx_get_stop_condition() == gmx_stop_cond_next))
1359 /* We need at least two global communication steps to pass
1360 * around the signal. We stop at a pair-list creation step
1361 * to allow for exact continuation, when possible.
1363 signals[eglsSTOPCOND].sig = 1;
1364 nsteps_stop = std::max(ir->nstlist, 2*nstSignalComm);
1366 else if (gmx_get_stop_condition() == gmx_stop_cond_next)
1368 /* Stop directly after the next global communication step.
1369 * This breaks exact continuation.
1371 signals[eglsSTOPCOND].sig = -1;
1372 nsteps_stop = nstSignalComm + 1;
1377 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1378 gmx_get_signal_name(), nsteps_stop);
1382 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1383 gmx_get_signal_name(), nsteps_stop);
1385 handled_stop_condition = static_cast<int>(gmx_get_stop_condition());
1387 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1388 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1389 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1391 /* Signal to terminate the run */
1392 signals[eglsSTOPCOND].sig = 1;
1395 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1397 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1400 if (bResetCountersHalfMaxH && MASTER(cr) &&
1401 elapsed_time > max_hours*60.0*60.0*0.495)
1403 /* Set flag that will communicate the signal to all ranks in the simulation */
1404 signals[eglsRESETCOUNTERS].sig = 1;
1407 /* In parallel we only have to check for checkpointing in steps
1408 * where we do global communication,
1409 * otherwise the other nodes don't know.
1411 const real cpt_period = mdrunOptions.checkpointOptions.period;
1412 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1415 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1416 signals[eglsCHKPT].set == 0)
1418 signals[eglsCHKPT].sig = 1;
1421 /* ######### START SECOND UPDATE STEP ################# */
1423 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1426 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1428 gmx_bool bIfRandomize;
1429 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1430 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1431 if (constr && bIfRandomize)
1433 constrain_velocities(step, nullptr,
1437 bCalcVir, do_log, do_ene);
1440 /* Box is changed in update() when we do pressure coupling,
1441 * but we should still use the old box for energy corrections and when
1442 * writing it to the energy file, so it matches the trajectory files for
1443 * the same timestep above. Make a copy in a separate array.
1445 copy_mat(state->box, lastbox);
1449 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1451 wallcycle_start(wcycle, ewcUPDATE);
1452 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1455 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1456 /* We can only do Berendsen coupling after we have summed
1457 * the kinetic energy or virial. Since the happens
1458 * in global_state after update, we should only do it at
1459 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1464 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1465 update_pcouple_before_coordinates(fplog, step, ir, state,
1466 parrinellorahmanMu, M,
1472 /* velocity half-step update */
1473 update_coords(step, ir, mdatoms, state, f, fcd,
1474 ekind, M, upd, etrtVELOCITY2,
1478 /* Above, initialize just copies ekinh into ekin,
1479 * it doesn't copy position (for VV),
1480 * and entire integrator for MD.
1483 if (ir->eI == eiVVAK)
1485 /* We probably only need md->homenr, not state->natoms */
1486 if (state->natoms > cbuf_nalloc)
1488 cbuf_nalloc = state->natoms;
1489 srenew(cbuf, cbuf_nalloc);
1491 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1494 update_coords(step, ir, mdatoms, state, f, fcd,
1495 ekind, M, upd, etrtPOSITION, cr, constr);
1496 wallcycle_stop(wcycle, ewcUPDATE);
1498 constrain_coordinates(step, &dvdl_constr, state,
1500 wcycle, upd, constr,
1501 bCalcVir, do_log, do_ene);
1502 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state,
1503 cr, nrnb, wcycle, upd, constr, do_log, do_ene);
1504 finish_update(ir, mdatoms,
1506 nrnb, wcycle, upd, constr);
1508 if (ir->eI == eiVVAK)
1510 /* erase F_EKIN and F_TEMP here? */
1511 /* just compute the kinetic energy at the half step to perform a trotter step */
1512 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1513 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1514 constr, &nullSignaller, lastbox,
1515 nullptr, &bSumEkinhOld,
1516 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1518 wallcycle_start(wcycle, ewcUPDATE);
1519 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1520 /* now we know the scaling, we can compute the positions again again */
1521 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1523 update_coords(step, ir, mdatoms, state, f, fcd,
1524 ekind, M, upd, etrtPOSITION, cr, constr);
1525 wallcycle_stop(wcycle, ewcUPDATE);
1527 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1528 /* are the small terms in the shake_vir here due
1529 * to numerical errors, or are they important
1530 * physically? I'm thinking they are just errors, but not completely sure.
1531 * For now, will call without actually constraining, constr=NULL*/
1532 finish_update(ir, mdatoms,
1534 nrnb, wcycle, upd, nullptr);
1538 /* this factor or 2 correction is necessary
1539 because half of the constraint force is removed
1540 in the vv step, so we have to double it. See
1541 the Redmine issue #1255. It is not yet clear
1542 if the factor of 2 is exact, or just a very
1543 good approximation, and this will be
1544 investigated. The next step is to see if this
1545 can be done adding a dhdl contribution from the
1546 rattle step, but this is somewhat more
1547 complicated with the current code. Will be
1548 investigated, hopefully for 4.6.3. However,
1549 this current solution is much better than
1550 having it completely wrong.
1552 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1556 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1561 /* Need to unshift here */
1562 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1565 if (vsite != nullptr)
1567 wallcycle_start(wcycle, ewcVSITECONSTR);
1568 if (graph != nullptr)
1570 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1572 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
1573 top->idef.iparams, top->idef.il,
1574 fr->ePBC, fr->bMolPBC, cr, state->box);
1576 if (graph != nullptr)
1578 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1580 wallcycle_stop(wcycle, ewcVSITECONSTR);
1583 /* ############## IF NOT VV, Calculate globals HERE ############ */
1584 /* With Leap-Frog we can skip compute_globals at
1585 * non-communication steps, but we need to calculate
1586 * the kinetic energy one step before communication.
1589 // Organize to do inter-simulation signalling on steps if
1590 // and when algorithms require it.
1591 bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1593 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1595 // Since we're already communicating at this step, we
1596 // can propagate intra-simulation signals. Note that
1597 // check_nstglobalcomm has the responsibility for
1598 // choosing the value of nstglobalcomm that is one way
1599 // bGStat becomes true, so we can't get into a
1600 // situation where e.g. checkpointing can't be
1602 bool doIntraSimSignal = true;
1603 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1605 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1606 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1609 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1610 (bGStat ? CGLO_GSTAT : 0)
1611 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1612 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1613 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1614 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1616 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1618 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1619 top_global, top, state,
1620 &shouldCheckNumberOfBondedInteractions);
1624 /* ############# END CALC EKIN AND PRESSURE ################# */
1626 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1627 the virial that should probably be addressed eventually. state->veta has better properies,
1628 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1629 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1631 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1633 /* Sum up the foreign energy and dhdl terms for md and sd.
1634 Currently done every step so that dhdl is correct in the .edr */
1635 sum_dhdl(enerd, state->lambda, ir->fepvals);
1638 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1639 pres, force_vir, shake_vir,
1643 /* ################# END UPDATE STEP 2 ################# */
1644 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1646 /* The coordinates (x) were unshifted in update */
1649 /* We will not sum ekinh_old,
1650 * so signal that we still have to do it.
1652 bSumEkinhOld = TRUE;
1657 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1659 /* use the directly determined last velocity, not actually the averaged half steps */
1660 if (bTrotter && ir->eI == eiVV)
1662 enerd->term[F_EKIN] = last_ekin;
1664 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1666 if (integratorHasConservedEnergyQuantity(ir))
1670 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1674 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1677 /* ######### END PREPARING EDR OUTPUT ########### */
1683 if (fplog && do_log && bDoExpanded)
1685 /* only needed if doing expanded ensemble */
1686 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1687 state_global->dfhist, state->fep_state, ir->nstlog, step);
1691 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1692 t, mdatoms->tmass, enerd, state,
1693 ir->fepvals, ir->expandedvals, lastbox,
1694 shake_vir, force_vir, total_vir, pres,
1695 ekind, mu_tot, constr);
1699 upd_mdebin_step(mdebin);
1702 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1703 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1705 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1707 eprNORMAL, mdebin, fcd, groups, &(ir->opts), awh.get());
1711 pull_print_output(ir->pull_work, step, t);
1714 if (do_per_step(step, ir->nstlog))
1716 if (fflush(fplog) != 0)
1718 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1724 /* Have to do this part _after_ outputting the logfile and the edr file */
1725 /* Gets written into the state at the beginning of next loop*/
1726 state->fep_state = lamnew;
1728 /* Print the remaining wall clock time for the run */
1729 if (isMasterSimMasterRank(ms, cr) &&
1730 (do_verbose || gmx_got_usr_signal()) &&
1735 fprintf(stderr, "\n");
1737 print_time(stderr, walltime_accounting, step, ir, cr);
1740 /* Ion/water position swapping.
1741 * Not done in last step since trajectory writing happens before this call
1742 * in the MD loop and exchanges would be lost anyway. */
1743 bNeedRepartition = FALSE;
1744 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1745 do_per_step(step, ir->swap->nstswap))
1747 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1748 bRerunMD ? rerun_fr.x : as_rvec_array(state->x.data()),
1749 bRerunMD ? rerun_fr.box : state->box,
1750 MASTER(cr) && mdrunOptions.verbose,
1753 if (bNeedRepartition && DOMAINDECOMP(cr))
1755 dd_collect_state(cr->dd, state, state_global);
1759 /* Replica exchange */
1763 bExchanged = replica_exchange(fplog, cr, ms, repl_ex,
1764 state_global, enerd,
1768 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1770 dd_partition_system(fplog, step, cr, TRUE, 1,
1771 state_global, top_global, ir,
1772 state, &f, mdAtoms, top, fr,
1774 nrnb, wcycle, FALSE);
1775 shouldCheckNumberOfBondedInteractions = true;
1776 update_realloc(upd, state->natoms);
1781 startingFromCheckpoint = false;
1783 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1784 /* With all integrators, except VV, we need to retain the pressure
1785 * at the current step for coupling at the next step.
1787 if ((state->flags & (1<<estPRES_PREV)) &&
1789 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1791 /* Store the pressure in t_state for pressure coupling
1792 * at the next MD step.
1794 copy_mat(pres, state->pres_prev);
1797 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1799 if ( (membed != nullptr) && (!bLastStep) )
1801 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1808 /* read next frame from input trajectory */
1809 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1814 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1818 cycles = wallcycle_stop(wcycle, ewcSTEP);
1819 if (DOMAINDECOMP(cr) && wcycle)
1821 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1824 if (!bRerunMD || !rerun_fr.bStep)
1826 /* increase the MD step number */
1831 /* TODO make a counter-reset module */
1832 /* If it is time to reset counters, set a flag that remains
1833 true until counters actually get reset */
1834 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1835 signals[eglsRESETCOUNTERS].set != 0)
1837 if (pme_loadbal_is_active(pme_loadbal))
1839 /* Do not permit counter reset while PME load
1840 * balancing is active. The only purpose for resetting
1841 * counters is to measure reliable performance data,
1842 * and that can't be done before balancing
1845 * TODO consider fixing this by delaying the reset
1846 * until after load balancing completes,
1847 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1848 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1849 "reset mdrun counters at step %" PRId64 ". Try "
1850 "resetting counters later in the run, e.g. with gmx "
1851 "mdrun -resetstep.", step);
1853 reset_all_counters(fplog, mdlog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1854 use_GPU(fr->nbv) ? fr->nbv : nullptr, fr->pmedata);
1855 wcycle_set_reset_counters(wcycle, -1);
1856 if (!thisRankHasDuty(cr, DUTY_PME))
1858 /* Tell our PME node to reset its counters */
1859 gmx_pme_send_resetcounters(cr, step);
1861 /* Correct max_hours for the elapsed time */
1862 max_hours -= elapsed_time/(60.0*60.0);
1863 /* If mdrun -maxh -resethway was active, it can only trigger once */
1864 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1865 /* Reset can only happen once, so clear the triggering flag. */
1866 signals[eglsRESETCOUNTERS].set = 0;
1869 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1870 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1873 /* End of main MD loop */
1875 /* Closing TNG files can include compressing data. Therefore it is good to do that
1876 * before stopping the time measurements. */
1877 mdoutf_tng_close(outf);
1879 /* Stop measuring walltime */
1880 walltime_accounting_end(walltime_accounting);
1882 if (bRerunMD && MASTER(cr))
1887 if (!thisRankHasDuty(cr, DUTY_PME))
1889 /* Tell the PME only node to finish */
1890 gmx_pme_send_finish(cr);
1895 if (ir->nstcalcenergy > 0 && !bRerunMD)
1897 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1898 eprAVER, mdebin, fcd, groups, &(ir->opts), awh.get());
1901 done_mdebin(mdebin);
1906 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
1909 done_shellfc(fplog, shellfc, step_rel);
1911 if (useReplicaExchange && MASTER(cr))
1913 print_replica_exchange_statistics(fplog, repl_ex);
1916 // Clean up swapcoords
1917 if (ir->eSwapCoords != eswapNO)
1919 finish_swapcoords(ir->swap);
1922 /* IMD cleanup, if bIMD is TRUE. */
1923 IMD_finalize(ir->bIMD, ir->imd);
1925 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1926 if (step_rel >= wcycle_get_reset_counters(wcycle) &&
1927 signals[eglsRESETCOUNTERS].set == 0 &&
1928 !bResetCountersHalfMaxH)
1930 walltime_accounting_set_valid_finish(walltime_accounting);
1933 destroy_enerdata(enerd);