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/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"
135 #include "replicaexchange.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_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 wallcycle_start(wcycle, ewcRUN);
181 walltime_accounting_reset_time(walltime_accounting, step);
182 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
185 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
187 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
188 * \param[in,out] globalState The global state container
189 * \param[in] constructVsites When true, vsite coordinates are constructed
190 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
191 * \param[in] idef Topology parameters, used for constructing vsites
192 * \param[in] timeStep Time step, used for constructing vsites
193 * \param[in] forceRec Force record, used for constructing vsites
194 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
195 * \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
197 static void prepareRerunState(const t_trxframe &rerunFrame,
198 t_state *globalState,
199 bool constructVsites,
200 const gmx_vsite_t *vsite,
203 const t_forcerec &forceRec,
205 gmx_bool *warnWhenNoV)
207 for (int i = 0; i < globalState->natoms; i++)
209 copy_rvec(rerunFrame.x[i], globalState->x[i]);
213 for (int i = 0; i < globalState->natoms; i++)
215 copy_rvec(rerunFrame.v[i], globalState->v[i]);
220 for (int i = 0; i < globalState->natoms; i++)
222 clear_rvec(globalState->v[i]);
226 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
227 " Ekin, temperature and pressure are incorrect,\n"
228 " the virial will be incorrect when constraints are present.\n"
230 *warnWhenNoV = FALSE;
233 copy_mat(rerunFrame.box, globalState->box);
237 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
241 /* Following is necessary because the graph may get out of sync
242 * with the coordinates if we only have every N'th coordinate set
244 mk_mshift(nullptr, graph, forceRec.ePBC, globalState->box, as_rvec_array(globalState->x.data()));
245 shift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
247 construct_vsites(vsite, as_rvec_array(globalState->x.data()), timeStep, as_rvec_array(globalState->v.data()),
248 idef.iparams, idef.il,
249 forceRec.ePBC, forceRec.bMolPBC, nullptr, globalState->box);
252 unshift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
257 void gmx::Integrator::do_md()
259 // TODO Historically, the EM and MD "integrators" used different
260 // names for the t_inputrec *parameter, but these must have the
261 // same name, now that it's a member of a struct. We use this ir
262 // alias to avoid a large ripple of nearly useless changes.
263 // t_inputrec is being replaced by IMdpOptionsProvider, so this
264 // will go away eventually.
265 t_inputrec *ir = inputrec;
266 gmx_mdoutf *outf = nullptr;
267 int64_t step, step_rel;
268 double t, t0, lam0[efptNR];
269 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
270 gmx_bool bNS, bNStList, bSimAnn, bStopCM,
271 bFirstStep, bInitStep, bLastStep = FALSE;
272 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
273 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
274 bForceUpdate = FALSE, bCPT;
275 gmx_bool bMasterState;
276 int force_flags, cglo_flags;
277 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
282 matrix parrinellorahmanMu, M;
284 gmx_repl_ex_t repl_ex = nullptr;
287 t_mdebin *mdebin = nullptr;
288 gmx_enerdata_t *enerd;
289 PaddedRVecVector f {};
290 gmx_global_stat_t gstat;
291 gmx_update_t *upd = nullptr;
292 t_graph *graph = nullptr;
293 gmx_groups_t *groups;
294 gmx_ekindata_t *ekind;
295 gmx_shellfc_t *shellfc;
296 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
297 gmx_bool bResetCountersHalfMaxH = FALSE;
298 gmx_bool bTemp, bPres, bTrotter;
300 rvec *cbuf = nullptr;
307 real saved_conserved_quantity = 0;
311 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
312 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
315 /* PME load balancing data for GPU kernels */
316 pme_load_balancing_t *pme_loadbal = nullptr;
317 gmx_bool bPMETune = FALSE;
318 gmx_bool bPMETunePrinting = FALSE;
321 gmx_bool bIMDstep = FALSE;
324 /* Temporary addition for FAHCORE checkpointing */
327 /* Domain decomposition could incorrectly miss a bonded
328 interaction, but checking for that requires a global
329 communication stage, which does not otherwise happen in DD
330 code. So we do that alongside the first global energy reduction
331 after a new DD is made. These variables handle whether the
332 check happens, and the result it returns. */
333 bool shouldCheckNumberOfBondedInteractions = false;
334 int totalNumberOfBondedInteractions = -1;
336 SimulationSignals signals;
337 // Most global communnication stages don't propagate mdrun
338 // signals, and will use this object to achieve that.
339 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
341 if (!mdrunOptions.writeConfout)
343 // This is on by default, and the main known use case for
344 // turning it off is for convenience in benchmarking, which is
345 // something that should not show up in the general user
347 GMX_LOG(mdlog.info).asParagraph().
348 appendText("The -noconfout functionality is deprecated, and may be removed in a future version.");
350 if (mdrunOptions.timingOptions.resetHalfway)
352 GMX_LOG(mdlog.info).asParagraph().
353 appendText("The -resethway functionality is deprecated, and may be removed in a future version.");
356 /* Signal to reset the counters half the simulation steps. */
357 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
359 /* Signal to reset the counters halfway the simulation time. */
360 bResetCountersHalfMaxH = (mdrunOptions.maximumHoursToRun > 0);
363 /* md-vv uses averaged full step velocities for T-control
364 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
365 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
366 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
368 const gmx_bool bRerunMD = mdrunOptions.rerun;
369 int nstglobalcomm = mdrunOptions.globalCommunicationInterval;
373 /* Since we don't know if the frames read are related in any way,
374 * rebuild the neighborlist at every step.
377 ir->nstcalcenergy = 1;
381 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
382 bGStatEveryStep = (nstglobalcomm == 1);
386 ir->nstxout_compressed = 0;
388 groups = &top_global->groups;
390 std::unique_ptr<EssentialDynamics> ed = nullptr;
391 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
393 /* Initialize essential dynamics sampling */
394 ed = init_edsam(opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm),
397 state_global, observablesHistory,
398 oenv, mdrunOptions.continuationOptions.appendFiles);
402 init_md(fplog, cr, outputProvider, ir, oenv, mdrunOptions,
403 &t, &t0, state_global, lam0,
404 nrnb, top_global, &upd, deform,
405 nfile, fnm, &outf, &mdebin,
406 force_vir, shake_vir, total_vir, pres, mu_tot, &bSimAnn, &vcm, wcycle);
408 /* Energy terms and groups */
410 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
413 /* Kinetic energy data */
415 init_ekindata(fplog, top_global, &(ir->opts), ekind);
416 /* Copy the cos acceleration to the groups struct */
417 ekind->cosacc.cos_accel = ir->cos_accel;
419 gstat = global_stat_init(ir);
421 /* Check for polarizable models and flexible constraints */
422 shellfc = init_shell_flexcon(fplog,
423 top_global, constr ? constr->numFlexibleConstraints() : 0,
424 ir->nstcalcenergy, DOMAINDECOMP(cr));
427 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
428 if ((io > 2000) && MASTER(cr))
431 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
436 /* Set up interactive MD (IMD) */
437 init_IMD(ir, cr, ms, top_global, fplog, ir->nstcalcenergy,
438 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
439 nfile, fnm, oenv, mdrunOptions);
441 // Local state only becomes valid now.
442 std::unique_ptr<t_state> stateInstance;
445 if (DOMAINDECOMP(cr))
447 top = dd_init_local_top(top_global);
449 stateInstance = compat::make_unique<t_state>();
450 state = stateInstance.get();
451 dd_init_local_state(cr->dd, state_global, state);
453 /* Distribute the charge groups over the nodes from the master node */
454 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
455 state_global, top_global, ir,
456 state, &f, mdAtoms, top, fr,
458 nrnb, nullptr, FALSE);
459 shouldCheckNumberOfBondedInteractions = true;
460 update_realloc(upd, state->natoms);
464 state_change_natoms(state_global, state_global->natoms);
465 /* We need to allocate one element extra, since we might use
466 * (unaligned) 4-wide SIMD loads to access rvec entries.
468 f.resize(gmx::paddedRVecVectorSize(state_global->natoms));
469 /* Copy the pointer to the global state */
470 state = state_global;
473 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
474 &graph, mdAtoms, constr, vsite, shellfc);
476 update_realloc(upd, state->natoms);
479 auto mdatoms = mdAtoms->mdatoms();
481 // NOTE: The global state is no longer used at this point.
482 // But state_global is still used as temporary storage space for writing
483 // the global state to file and potentially for replica exchange.
484 // (Global topology should persist.)
486 update_mdatoms(mdatoms, state->lambda[efptMASS]);
488 const ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
489 bool startingFromCheckpoint = continuationOptions.startedFromCheckpoint;
493 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
498 if (startingFromCheckpoint)
500 /* Update mdebin with energy history if appending to output files */
501 if (continuationOptions.appendFiles)
503 restore_energyhistory_from_state(mdebin, observablesHistory->energyHistory.get());
505 else if (observablesHistory->energyHistory != nullptr)
507 /* We might have read an energy history from checkpoint.
508 * As we are not appending, we want to restart the statistics.
509 * Free the allocated memory and reset the counts.
511 observablesHistory->energyHistory = {};
514 if (observablesHistory->energyHistory == nullptr)
516 observablesHistory->energyHistory = compat::make_unique<energyhistory_t>();
518 /* Set the initial energy history in state by updating once */
519 update_energyhistory(observablesHistory->energyHistory.get(), mdebin);
522 // TODO: Remove this by converting AWH into a ForceProvider
523 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms, startingFromCheckpoint,
525 opt2fn("-awh", nfile, fnm), ir->pull_work);
527 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
528 if (useReplicaExchange && MASTER(cr))
530 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir,
533 /* PME tuning is only supported in the Verlet scheme, with PME for
534 * Coulomb. It is not supported with only LJ PME, or for
536 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !bRerunMD &&
537 !mdrunOptions.reproducible && ir->cutoff_scheme != ecutsGROUP);
540 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
541 fr->ic, fr->nbv->listParams.get(), fr->pmedata, use_GPU(fr->nbv),
545 if (!ir->bContinuation && !bRerunMD)
547 if (state->flags & (1 << estV))
549 /* Set the velocities of vsites, shells and frozen atoms to zero */
550 for (i = 0; i < mdatoms->homenr; i++)
552 if (mdatoms->ptype[i] == eptVSite ||
553 mdatoms->ptype[i] == eptShell)
555 clear_rvec(state->v[i]);
557 else if (mdatoms->cFREEZE)
559 for (m = 0; m < DIM; m++)
561 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
572 /* Constrain the initial coordinates and velocities */
573 do_constrain_first(fplog, constr, ir, mdatoms, state);
577 /* Construct the virtual sites for the initial configuration */
578 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, nullptr,
579 top->idef.iparams, top->idef.il,
580 fr->ePBC, fr->bMolPBC, cr, state->box);
584 if (ir->efep != efepNO)
586 /* Set free energy calculation frequency as the greatest common
587 * denominator of nstdhdl and repl_ex_nst. */
588 nstfep = ir->fepvals->nstdhdl;
591 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
593 if (useReplicaExchange)
595 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
599 /* Be REALLY careful about what flags you set here. You CANNOT assume
600 * this is the first step, since we might be restarting from a checkpoint,
601 * and in that case we should not do any modifications to the state.
603 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
605 if (continuationOptions.haveReadEkin)
607 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
610 cglo_flags = (CGLO_INITIALIZATION | CGLO_TEMPERATURE | CGLO_GSTAT
611 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
612 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
613 | (continuationOptions.haveReadEkin ? CGLO_READEKIN : 0));
615 bSumEkinhOld = FALSE;
616 /* To minimize communication, compute_globals computes the COM velocity
617 * and the kinetic energy for the velocities without COM motion removed.
618 * Thus to get the kinetic energy without the COM contribution, we need
619 * to call compute_globals twice.
621 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
623 int cglo_flags_iteration = cglo_flags;
624 if (bStopCM && cgloIteration == 0)
626 cglo_flags_iteration |= CGLO_STOPCM;
627 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
629 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
630 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
631 constr, &nullSignaller, state->box,
632 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags_iteration
633 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
635 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
636 top_global, top, state,
637 &shouldCheckNumberOfBondedInteractions);
638 if (ir->eI == eiVVAK)
640 /* a second call to get the half step temperature initialized as well */
641 /* we do the same call as above, but turn the pressure off -- internally to
642 compute_globals, this is recognized as a velocity verlet half-step
643 kinetic energy calculation. This minimized excess variables, but
644 perhaps loses some logic?*/
646 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
647 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
648 constr, &nullSignaller, state->box,
649 nullptr, &bSumEkinhOld,
650 cglo_flags & ~CGLO_PRESSURE);
653 /* Calculate the initial half step temperature, and save the ekinh_old */
654 if (!continuationOptions.startedFromCheckpoint)
656 for (i = 0; (i < ir->opts.ngtc); i++)
658 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
662 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
663 temperature control */
664 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
668 if (!ir->bContinuation)
670 if (constr && ir->eConstrAlg == econtLINCS)
673 "RMS relative constraint deviation after constraining: %.2e\n",
676 if (EI_STATE_VELOCITY(ir->eI))
678 real temp = enerd->term[F_TEMP];
681 /* Result of Ekin averaged over velocities of -half
682 * and +half step, while we only have -half step here.
686 fprintf(fplog, "Initial temperature: %g K\n", temp);
692 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
693 " input trajectory '%s'\n\n",
694 *(top_global->name), opt2fn("-rerun", nfile, fnm));
695 if (mdrunOptions.verbose)
697 fprintf(stderr, "Calculated time to finish depends on nsteps from "
698 "run input file,\nwhich may not correspond to the time "
699 "needed to process input trajectory.\n\n");
705 fprintf(stderr, "starting mdrun '%s'\n",
706 *(top_global->name));
709 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
713 sprintf(tbuf, "%s", "infinite");
715 if (ir->init_step > 0)
717 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
718 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
719 gmx_step_str(ir->init_step, sbuf2),
720 ir->init_step*ir->delta_t);
724 fprintf(stderr, "%s steps, %s ps.\n",
725 gmx_step_str(ir->nsteps, sbuf), tbuf);
728 fprintf(fplog, "\n");
731 walltime_accounting_start_time(walltime_accounting);
732 wallcycle_start(wcycle, ewcRUN);
733 print_start(fplog, cr, walltime_accounting, "mdrun");
735 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
737 chkpt_ret = fcCheckPointParallel( cr->nodeid,
741 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
745 /***********************************************************
749 ************************************************************/
751 /* if rerunMD then read coordinates and velocities from input trajectory */
754 if (getenv("GMX_FORCE_UPDATE"))
762 bLastStep = !read_first_frame(oenv, &status,
763 opt2fn("-rerun", nfile, fnm),
764 &rerun_fr, TRX_NEED_X | TRX_READ_V);
765 if (rerun_fr.natoms != top_global->natoms)
768 "Number of atoms in trajectory (%d) does not match the "
769 "run input file (%d)\n",
770 rerun_fr.natoms, top_global->natoms);
772 if (ir->ePBC != epbcNONE)
776 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f does not contain a box, while pbc is used", rerun_fr.step, rerun_fr.time);
778 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
780 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
787 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
790 if (ir->ePBC != epbcNONE)
792 /* Set the shift vectors.
793 * Necessary here when have a static box different from the tpr box.
795 calc_shifts(rerun_fr.box, fr->shift_vec);
800 /* Skip the first Nose-Hoover integration when we get the state from tpx */
801 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
802 bSumEkinhOld = FALSE;
804 bNeedRepartition = FALSE;
806 bool simulationsShareState = false;
807 int nstSignalComm = nstglobalcomm;
809 // TODO This implementation of ensemble orientation restraints is nasty because
810 // a user can't just do multi-sim with single-sim orientation restraints.
811 bool usingEnsembleRestraints = (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
812 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
814 // Replica exchange, ensemble restraints and AWH need all
815 // simulations to remain synchronized, so they need
816 // checkpoints and stop conditions to act on the same step, so
817 // the propagation of such signals must take place between
818 // simulations, not just within simulations.
819 // TODO: Make algorithm initializers set these flags.
820 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
821 bool resetCountersIsLocal = true;
822 signals[eglsCHKPT] = SimulationSignal(!simulationsShareState);
823 signals[eglsSTOPCOND] = SimulationSignal(!simulationsShareState);
824 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
826 if (simulationsShareState)
828 // Inter-simulation signal communication does not need to happen
829 // often, so we use a minimum of 200 steps to reduce overhead.
830 const int c_minimumInterSimulationSignallingInterval = 200;
831 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1)/nstglobalcomm)*nstglobalcomm;
835 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion = (DOMAINDECOMP(cr) ? DdOpenBalanceRegionBeforeForceComputation::yes : DdOpenBalanceRegionBeforeForceComputation::no);
836 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion = (DOMAINDECOMP(cr) ? DdCloseBalanceRegionAfterForceComputation::yes : DdCloseBalanceRegionAfterForceComputation::no);
838 step = ir->init_step;
841 // TODO extract this to new multi-simulation module
842 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
844 if (!multisim_int_all_are_equal(ms, ir->nsteps))
846 GMX_LOG(mdlog.warning).appendText(
847 "Note: The number of steps is not consistent across multi simulations,\n"
848 "but we are proceeding anyway!");
850 if (!multisim_int_all_are_equal(ms, ir->init_step))
852 GMX_LOG(mdlog.warning).appendText(
853 "Note: The initial step is not consistent across multi simulations,\n"
854 "but we are proceeding anyway!");
858 /* and stop now if we should */
859 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
863 /* Determine if this is a neighbor search step */
864 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
866 if (bPMETune && bNStList)
868 /* PME grid + cut-off optimization with GPUs or PME nodes */
869 pme_loadbal_do(pme_loadbal, cr,
870 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
878 wallcycle_start(wcycle, ewcSTEP);
884 step = rerun_fr.step;
885 step_rel = step - ir->init_step;
898 bLastStep = (step_rel == ir->nsteps);
899 t = t0 + step*ir->delta_t;
902 // TODO Refactor this, so that nstfep does not need a default value of zero
903 if (ir->efep != efepNO || ir->bSimTemp)
905 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
906 requiring different logic. */
911 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
916 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
918 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
919 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
920 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
921 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
924 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep &&
925 do_per_step(step, replExParams.exchangeInterval));
929 update_annealing_target_temp(ir, t, upd);
932 if (bRerunMD && MASTER(cr))
934 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
935 if (constructVsites && DOMAINDECOMP(cr))
937 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
939 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top->idef, ir->delta_t, *fr, graph, &bRerunWarnNoV);
942 /* Stop Center of Mass motion */
943 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
947 /* for rerun MD always do Neighbour Searching */
948 bNS = (bFirstStep || ir->nstlist != 0);
952 /* Determine whether or not to do Neighbour Searching */
953 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
956 /* < 0 means stop at next step, > 0 means stop at next NS step */
957 if ( (signals[eglsSTOPCOND].set < 0) ||
958 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
963 /* do_log triggers energy and virial calculation. Because this leads
964 * to different code paths, forces can be different. Thus for exact
965 * continuation we should avoid extra log output.
966 * Note that the || bLastStep can result in non-exact continuation
967 * beyond the last step. But we don't consider that to be an issue.
969 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
970 do_verbose = mdrunOptions.verbose &&
971 (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep || bRerunMD);
973 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
981 bMasterState = FALSE;
982 /* Correct the new box if it is too skewed */
983 if (inputrecDynamicBox(ir))
985 if (correct_box(fplog, step, state->box, graph))
990 if (DOMAINDECOMP(cr) && bMasterState)
992 dd_collect_state(cr->dd, state, state_global);
996 if (DOMAINDECOMP(cr))
998 /* Repartition the domain decomposition */
999 dd_partition_system(fplog, mdlog, step, cr,
1000 bMasterState, nstglobalcomm,
1001 state_global, top_global, ir,
1002 state, &f, mdAtoms, top, fr,
1005 do_verbose && !bPMETunePrinting);
1006 shouldCheckNumberOfBondedInteractions = true;
1007 update_realloc(upd, state->natoms);
1011 if (MASTER(cr) && do_log)
1013 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
1016 if (ir->efep != efepNO)
1018 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1021 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1024 /* We need the kinetic energy at minus the half step for determining
1025 * the full step kinetic energy and possibly for T-coupling.*/
1026 /* This may not be quite working correctly yet . . . . */
1027 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1028 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1029 constr, &nullSignaller, state->box,
1030 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1031 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1032 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1033 top_global, top, state,
1034 &shouldCheckNumberOfBondedInteractions);
1036 clear_mat(force_vir);
1038 /* We write a checkpoint at this MD step when:
1039 * either at an NS step when we signalled through gs,
1040 * or at the last step (but not when we do not want confout),
1041 * but never at the first step or with rerun.
1043 bCPT = ((((signals[eglsCHKPT].set != 0) && (bNS || ir->nstlist == 0)) ||
1044 (bLastStep && mdrunOptions.writeConfout)) &&
1045 step > ir->init_step && !bRerunMD);
1048 signals[eglsCHKPT].set = 0;
1051 /* Determine the energy and pressure:
1052 * at nstcalcenergy steps and at energy output steps (set below).
1054 if (EI_VV(ir->eI) && (!bInitStep))
1056 /* for vv, the first half of the integration actually corresponds
1057 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1058 but the virial needs to be calculated on both the current step and the 'next' step. Future
1059 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1061 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1062 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1063 bCalcVir = bCalcEnerStep ||
1064 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1068 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1069 bCalcVir = bCalcEnerStep ||
1070 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1072 bCalcEner = bCalcEnerStep;
1074 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1076 if (do_ene || do_log || bDoReplEx)
1082 /* Do we need global communication ? */
1083 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1084 do_per_step(step, nstglobalcomm) ||
1085 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1087 force_flags = (GMX_FORCE_STATECHANGED |
1088 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1089 GMX_FORCE_ALLFORCES |
1090 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1091 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1092 (bDoFEP ? GMX_FORCE_DHDL : 0)
1097 /* Now is the time to relax the shells */
1098 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
1099 enforcedRotation, step,
1100 ir, bNS, force_flags, top,
1102 state, f, force_vir, mdatoms,
1103 nrnb, wcycle, graph, groups,
1104 shellfc, fr, t, mu_tot,
1106 ddOpenBalanceRegion, ddCloseBalanceRegion);
1110 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
1111 (or the AWH update will be performed twice for one step when continuing). It would be best to
1112 call this update function from do_md_trajectory_writing but that would occur after do_force.
1113 One would have to divide the update_awh function into one function applying the AWH force
1114 and one doing the AWH bias update. The update AWH bias function could then be called after
1115 do_md_trajectory_writing (then containing update_awh_history).
1116 The checkpointing will in the future probably moved to the start of the md loop which will
1117 rid of this issue. */
1118 if (awh && bCPT && MASTER(cr))
1120 awh->updateHistory(state_global->awhHistory.get());
1123 /* The coordinates (x) are shifted (to get whole molecules)
1125 * This is parallellized as well, and does communication too.
1126 * Check comments in sim_util.c
1128 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation,
1129 step, nrnb, wcycle, top, groups,
1130 state->box, state->x, &state->hist,
1131 f, force_vir, mdatoms, enerd, fcd,
1132 state->lambda, graph,
1133 fr, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
1134 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1135 ddOpenBalanceRegion, ddCloseBalanceRegion);
1138 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1139 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1141 rvec *vbuf = nullptr;
1143 wallcycle_start(wcycle, ewcUPDATE);
1144 if (ir->eI == eiVV && bInitStep)
1146 /* if using velocity verlet with full time step Ekin,
1147 * take the first half step only to compute the
1148 * virial for the first step. From there,
1149 * revert back to the initial coordinates
1150 * so that the input is actually the initial step.
1152 snew(vbuf, state->natoms);
1153 copy_rvecn(as_rvec_array(state->v.data()), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1157 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1158 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1161 update_coords(step, ir, mdatoms, state, f, fcd,
1162 ekind, M, upd, etrtVELOCITY1,
1165 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1167 wallcycle_stop(wcycle, ewcUPDATE);
1168 constrain_velocities(step, nullptr,
1172 bCalcVir, do_log, do_ene);
1173 wallcycle_start(wcycle, ewcUPDATE);
1177 /* Need to unshift here if a do_force has been
1178 called in the previous step */
1179 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1181 /* if VV, compute the pressure and constraints */
1182 /* For VV2, we strictly only need this if using pressure
1183 * control, but we really would like to have accurate pressures
1185 * Think about ways around this in the future?
1186 * For now, keep this choice in comments.
1188 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1189 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1191 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1192 if (bCalcEner && ir->eI == eiVVAK)
1194 bSumEkinhOld = TRUE;
1196 /* for vv, the first half of the integration actually corresponds to the previous step.
1197 So we need information from the last step in the first half of the integration */
1198 if (bGStat || do_per_step(step-1, nstglobalcomm))
1200 wallcycle_stop(wcycle, ewcUPDATE);
1201 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1202 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1203 constr, &nullSignaller, state->box,
1204 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1205 (bGStat ? CGLO_GSTAT : 0)
1207 | (bTemp ? CGLO_TEMPERATURE : 0)
1208 | (bPres ? CGLO_PRESSURE : 0)
1209 | (bPres ? CGLO_CONSTRAINT : 0)
1210 | (bStopCM ? CGLO_STOPCM : 0)
1211 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1214 /* explanation of above:
1215 a) We compute Ekin at the full time step
1216 if 1) we are using the AveVel Ekin, and it's not the
1217 initial step, or 2) if we are using AveEkin, but need the full
1218 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1219 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1220 EkinAveVel because it's needed for the pressure */
1221 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1222 top_global, top, state,
1223 &shouldCheckNumberOfBondedInteractions);
1224 wallcycle_start(wcycle, ewcUPDATE);
1226 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1231 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1232 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1234 /* TODO This is only needed when we're about to write
1235 * a checkpoint, because we use it after the restart
1236 * (in a kludge?). But what should we be doing if
1237 * startingFromCheckpoint or bInitStep are true? */
1238 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1240 copy_mat(shake_vir, state->svir_prev);
1241 copy_mat(force_vir, state->fvir_prev);
1243 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1245 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1246 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1247 enerd->term[F_EKIN] = trace(ekind->ekin);
1250 else if (bExchanged)
1252 wallcycle_stop(wcycle, ewcUPDATE);
1253 /* We need the kinetic energy at minus the half step for determining
1254 * the full step kinetic energy and possibly for T-coupling.*/
1255 /* This may not be quite working correctly yet . . . . */
1256 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1257 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1258 constr, &nullSignaller, state->box,
1259 nullptr, &bSumEkinhOld,
1260 CGLO_GSTAT | CGLO_TEMPERATURE);
1261 wallcycle_start(wcycle, ewcUPDATE);
1264 /* if it's the initial step, we performed this first step just to get the constraint virial */
1265 if (ir->eI == eiVV && bInitStep)
1267 copy_rvecn(vbuf, as_rvec_array(state->v.data()), 0, state->natoms);
1270 wallcycle_stop(wcycle, ewcUPDATE);
1273 /* compute the conserved quantity */
1276 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1279 last_ekin = enerd->term[F_EKIN];
1281 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1283 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1285 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1286 if (ir->efep != efepNO && !bRerunMD)
1288 sum_dhdl(enerd, state->lambda, ir->fepvals);
1292 /* ######## END FIRST UPDATE STEP ############## */
1293 /* ######## If doing VV, we now have v(dt) ###### */
1296 /* perform extended ensemble sampling in lambda - we don't
1297 actually move to the new state before outputting
1298 statistics, but if performing simulated tempering, we
1299 do update the velocities and the tau_t. */
1301 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, as_rvec_array(state->v.data()), mdatoms);
1302 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1305 copy_df_history(state_global->dfhist, state->dfhist);
1309 /* Now we have the energies and forces corresponding to the
1310 * coordinates at time t. We must output all of this before
1313 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1314 ir, state, state_global, observablesHistory,
1316 outf, mdebin, ekind, f,
1318 bCPT, bRerunMD, bLastStep,
1319 mdrunOptions.writeConfout,
1321 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1322 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, as_rvec_array(state->x.data()), ir, t, wcycle);
1324 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1325 if (startingFromCheckpoint && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1327 copy_mat(state->svir_prev, shake_vir);
1328 copy_mat(state->fvir_prev, force_vir);
1331 double secondsSinceStart = walltime_accounting_get_time_since_start(walltime_accounting);
1333 /* Check whether everything is still allright */
1334 if ((static_cast<int>(gmx_get_stop_condition()) > handled_stop_condition)
1340 int nsteps_stop = -1;
1342 /* this just makes signals[].sig compatible with the hack
1343 of sending signals around by MPI_Reduce together with
1345 if ((gmx_get_stop_condition() == gmx_stop_cond_next_ns) ||
1346 (mdrunOptions.reproducible &&
1347 gmx_get_stop_condition() == gmx_stop_cond_next))
1349 /* We need at least two global communication steps to pass
1350 * around the signal. We stop at a pair-list creation step
1351 * to allow for exact continuation, when possible.
1353 signals[eglsSTOPCOND].sig = 1;
1354 nsteps_stop = std::max(ir->nstlist, 2*nstSignalComm);
1356 else if (gmx_get_stop_condition() == gmx_stop_cond_next)
1358 /* Stop directly after the next global communication step.
1359 * This breaks exact continuation.
1361 signals[eglsSTOPCOND].sig = -1;
1362 nsteps_stop = nstSignalComm + 1;
1367 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1368 gmx_get_signal_name(), nsteps_stop);
1372 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1373 gmx_get_signal_name(), nsteps_stop);
1375 handled_stop_condition = static_cast<int>(gmx_get_stop_condition());
1377 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1378 (mdrunOptions.maximumHoursToRun > 0 &&
1379 secondsSinceStart > mdrunOptions.maximumHoursToRun*60.0*60.0*0.99) &&
1380 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1382 /* Signal to terminate the run */
1383 signals[eglsSTOPCOND].sig = 1;
1386 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",
1387 gmx_step_str(step, sbuf), mdrunOptions.maximumHoursToRun*0.99);
1389 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",
1390 gmx_step_str(step, sbuf), mdrunOptions.maximumHoursToRun*0.99);
1393 if (bResetCountersHalfMaxH && MASTER(cr) &&
1394 secondsSinceStart > mdrunOptions.maximumHoursToRun*60.0*60.0*0.495)
1396 /* Set flag that will communicate the signal to all ranks in the simulation */
1397 signals[eglsRESETCOUNTERS].sig = 1;
1400 /* In parallel we only have to check for checkpointing in steps
1401 * where we do global communication,
1402 * otherwise the other nodes don't know.
1404 const real cpt_period = mdrunOptions.checkpointOptions.period;
1405 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1408 secondsSinceStart >= nchkpt*cpt_period*60.0)) &&
1409 signals[eglsCHKPT].set == 0)
1411 signals[eglsCHKPT].sig = 1;
1414 /* ######### START SECOND UPDATE STEP ################# */
1416 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1419 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1421 gmx_bool bIfRandomize;
1422 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1423 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1424 if (constr && bIfRandomize)
1426 constrain_velocities(step, nullptr,
1430 bCalcVir, do_log, do_ene);
1433 /* Box is changed in update() when we do pressure coupling,
1434 * but we should still use the old box for energy corrections and when
1435 * writing it to the energy file, so it matches the trajectory files for
1436 * the same timestep above. Make a copy in a separate array.
1438 copy_mat(state->box, lastbox);
1442 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1444 wallcycle_start(wcycle, ewcUPDATE);
1445 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1448 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1449 /* We can only do Berendsen coupling after we have summed
1450 * the kinetic energy or virial. Since the happens
1451 * in global_state after update, we should only do it at
1452 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1457 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1458 update_pcouple_before_coordinates(fplog, step, ir, state,
1459 parrinellorahmanMu, M,
1465 /* velocity half-step update */
1466 update_coords(step, ir, mdatoms, state, f, fcd,
1467 ekind, M, upd, etrtVELOCITY2,
1471 /* Above, initialize just copies ekinh into ekin,
1472 * it doesn't copy position (for VV),
1473 * and entire integrator for MD.
1476 if (ir->eI == eiVVAK)
1478 /* We probably only need md->homenr, not state->natoms */
1479 if (state->natoms > cbuf_nalloc)
1481 cbuf_nalloc = state->natoms;
1482 srenew(cbuf, cbuf_nalloc);
1484 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1487 update_coords(step, ir, mdatoms, state, f, fcd,
1488 ekind, M, upd, etrtPOSITION, cr, constr);
1489 wallcycle_stop(wcycle, ewcUPDATE);
1491 constrain_coordinates(step, &dvdl_constr, state,
1493 wcycle, upd, constr,
1494 bCalcVir, do_log, do_ene);
1495 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state,
1496 cr, nrnb, wcycle, upd, constr, do_log, do_ene);
1497 finish_update(ir, mdatoms,
1499 nrnb, wcycle, upd, constr);
1501 if (ir->eI == eiVVAK)
1503 /* erase F_EKIN and F_TEMP here? */
1504 /* just compute the kinetic energy at the half step to perform a trotter step */
1505 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1506 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1507 constr, &nullSignaller, lastbox,
1508 nullptr, &bSumEkinhOld,
1509 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1511 wallcycle_start(wcycle, ewcUPDATE);
1512 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1513 /* now we know the scaling, we can compute the positions again again */
1514 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1516 update_coords(step, ir, mdatoms, state, f, fcd,
1517 ekind, M, upd, etrtPOSITION, cr, constr);
1518 wallcycle_stop(wcycle, ewcUPDATE);
1520 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1521 /* are the small terms in the shake_vir here due
1522 * to numerical errors, or are they important
1523 * physically? I'm thinking they are just errors, but not completely sure.
1524 * For now, will call without actually constraining, constr=NULL*/
1525 finish_update(ir, mdatoms,
1527 nrnb, wcycle, upd, nullptr);
1531 /* this factor or 2 correction is necessary
1532 because half of the constraint force is removed
1533 in the vv step, so we have to double it. See
1534 the Redmine issue #1255. It is not yet clear
1535 if the factor of 2 is exact, or just a very
1536 good approximation, and this will be
1537 investigated. The next step is to see if this
1538 can be done adding a dhdl contribution from the
1539 rattle step, but this is somewhat more
1540 complicated with the current code. Will be
1541 investigated, hopefully for 4.6.3. However,
1542 this current solution is much better than
1543 having it completely wrong.
1545 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1549 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1554 /* Need to unshift here */
1555 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1558 if (vsite != nullptr)
1560 wallcycle_start(wcycle, ewcVSITECONSTR);
1561 if (graph != nullptr)
1563 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1565 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
1566 top->idef.iparams, top->idef.il,
1567 fr->ePBC, fr->bMolPBC, cr, state->box);
1569 if (graph != nullptr)
1571 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1573 wallcycle_stop(wcycle, ewcVSITECONSTR);
1576 /* ############## IF NOT VV, Calculate globals HERE ############ */
1577 /* With Leap-Frog we can skip compute_globals at
1578 * non-communication steps, but we need to calculate
1579 * the kinetic energy one step before communication.
1582 // Organize to do inter-simulation signalling on steps if
1583 // and when algorithms require it.
1584 bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1586 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1588 // Since we're already communicating at this step, we
1589 // can propagate intra-simulation signals. Note that
1590 // check_nstglobalcomm has the responsibility for
1591 // choosing the value of nstglobalcomm that is one way
1592 // bGStat becomes true, so we can't get into a
1593 // situation where e.g. checkpointing can't be
1595 bool doIntraSimSignal = true;
1596 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1598 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1599 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1602 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1603 (bGStat ? CGLO_GSTAT : 0)
1604 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1605 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1606 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1607 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1609 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1611 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1612 top_global, top, state,
1613 &shouldCheckNumberOfBondedInteractions);
1617 /* ############# END CALC EKIN AND PRESSURE ################# */
1619 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1620 the virial that should probably be addressed eventually. state->veta has better properies,
1621 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1622 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1624 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1626 /* Sum up the foreign energy and dhdl terms for md and sd.
1627 Currently done every step so that dhdl is correct in the .edr */
1628 sum_dhdl(enerd, state->lambda, ir->fepvals);
1631 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1632 pres, force_vir, shake_vir,
1636 /* ################# END UPDATE STEP 2 ################# */
1637 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1639 /* The coordinates (x) were unshifted in update */
1642 /* We will not sum ekinh_old,
1643 * so signal that we still have to do it.
1645 bSumEkinhOld = TRUE;
1650 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1652 /* use the directly determined last velocity, not actually the averaged half steps */
1653 if (bTrotter && ir->eI == eiVV)
1655 enerd->term[F_EKIN] = last_ekin;
1657 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1659 if (integratorHasConservedEnergyQuantity(ir))
1663 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1667 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1670 /* ######### END PREPARING EDR OUTPUT ########### */
1676 if (fplog && do_log && bDoExpanded)
1678 /* only needed if doing expanded ensemble */
1679 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1680 state_global->dfhist, state->fep_state, ir->nstlog, step);
1684 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1685 t, mdatoms->tmass, enerd, state,
1686 ir->fepvals, ir->expandedvals, lastbox,
1687 shake_vir, force_vir, total_vir, pres,
1688 ekind, mu_tot, constr);
1692 upd_mdebin_step(mdebin);
1695 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1696 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1698 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1700 eprNORMAL, mdebin, fcd, groups, &(ir->opts), awh.get());
1704 pull_print_output(ir->pull_work, step, t);
1707 if (do_per_step(step, ir->nstlog))
1709 if (fflush(fplog) != 0)
1711 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1717 /* Have to do this part _after_ outputting the logfile and the edr file */
1718 /* Gets written into the state at the beginning of next loop*/
1719 state->fep_state = lamnew;
1721 /* Print the remaining wall clock time for the run */
1722 if (isMasterSimMasterRank(ms, cr) &&
1723 (do_verbose || gmx_got_usr_signal()) &&
1728 fprintf(stderr, "\n");
1730 print_time(stderr, walltime_accounting, step, ir, cr);
1733 /* Ion/water position swapping.
1734 * Not done in last step since trajectory writing happens before this call
1735 * in the MD loop and exchanges would be lost anyway. */
1736 bNeedRepartition = FALSE;
1737 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1738 do_per_step(step, ir->swap->nstswap))
1740 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1741 bRerunMD ? rerun_fr.x : as_rvec_array(state->x.data()),
1742 bRerunMD ? rerun_fr.box : state->box,
1743 MASTER(cr) && mdrunOptions.verbose,
1746 if (bNeedRepartition && DOMAINDECOMP(cr))
1748 dd_collect_state(cr->dd, state, state_global);
1752 /* Replica exchange */
1756 bExchanged = replica_exchange(fplog, cr, ms, repl_ex,
1757 state_global, enerd,
1761 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1763 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1,
1764 state_global, top_global, ir,
1765 state, &f, mdAtoms, top, fr,
1767 nrnb, wcycle, FALSE);
1768 shouldCheckNumberOfBondedInteractions = true;
1769 update_realloc(upd, state->natoms);
1774 startingFromCheckpoint = false;
1776 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1777 /* With all integrators, except VV, we need to retain the pressure
1778 * at the current step for coupling at the next step.
1780 if ((state->flags & (1<<estPRES_PREV)) &&
1782 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1784 /* Store the pressure in t_state for pressure coupling
1785 * at the next MD step.
1787 copy_mat(pres, state->pres_prev);
1790 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1792 if ( (membed != nullptr) && (!bLastStep) )
1794 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1801 /* read next frame from input trajectory */
1802 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1807 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1811 cycles = wallcycle_stop(wcycle, ewcSTEP);
1812 if (DOMAINDECOMP(cr) && wcycle)
1814 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1817 if (!bRerunMD || !rerun_fr.bStep)
1819 /* increase the MD step number */
1824 /* TODO make a counter-reset module */
1825 /* If it is time to reset counters, set a flag that remains
1826 true until counters actually get reset */
1827 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1828 signals[eglsRESETCOUNTERS].set != 0)
1830 if (pme_loadbal_is_active(pme_loadbal))
1832 /* Do not permit counter reset while PME load
1833 * balancing is active. The only purpose for resetting
1834 * counters is to measure reliable performance data,
1835 * and that can't be done before balancing
1838 * TODO consider fixing this by delaying the reset
1839 * until after load balancing completes,
1840 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1841 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1842 "reset mdrun counters at step %" PRId64 ". Try "
1843 "resetting counters later in the run, e.g. with gmx "
1844 "mdrun -resetstep.", step);
1846 reset_all_counters(fplog, mdlog, cr, step, wcycle, nrnb, walltime_accounting,
1847 use_GPU(fr->nbv) ? fr->nbv : nullptr, fr->pmedata);
1848 wcycle_set_reset_counters(wcycle, -1);
1849 if (!thisRankHasDuty(cr, DUTY_PME))
1851 /* Tell our PME node to reset its counters */
1852 gmx_pme_send_resetcounters(cr, step);
1854 /* If mdrun -maxh -resethway was active, it can only trigger once */
1855 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1856 /* Reset can only happen once, so clear the triggering flag. */
1857 signals[eglsRESETCOUNTERS].set = 0;
1860 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1861 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1864 /* End of main MD loop */
1866 /* Closing TNG files can include compressing data. Therefore it is good to do that
1867 * before stopping the time measurements. */
1868 mdoutf_tng_close(outf);
1870 /* Stop measuring walltime */
1871 walltime_accounting_end_time(walltime_accounting);
1873 if (bRerunMD && MASTER(cr))
1878 if (!thisRankHasDuty(cr, DUTY_PME))
1880 /* Tell the PME only node to finish */
1881 gmx_pme_send_finish(cr);
1886 if (ir->nstcalcenergy > 0 && !bRerunMD)
1888 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1889 eprAVER, mdebin, fcd, groups, &(ir->opts), awh.get());
1892 done_mdebin(mdebin);
1897 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
1900 done_shellfc(fplog, shellfc, step_rel);
1902 if (useReplicaExchange && MASTER(cr))
1904 print_replica_exchange_statistics(fplog, repl_ex);
1907 // Clean up swapcoords
1908 if (ir->eSwapCoords != eswapNO)
1910 finish_swapcoords(ir->swap);
1913 /* IMD cleanup, if bIMD is TRUE. */
1914 IMD_finalize(ir->bIMD, ir->imd);
1916 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1917 if (step_rel >= wcycle_get_reset_counters(wcycle) &&
1918 signals[eglsRESETCOUNTERS].set == 0 &&
1919 !bResetCountersHalfMaxH)
1921 walltime_accounting_set_valid_finish(walltime_accounting);
1924 destroy_enerdata(enerd);