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 "thread_mpi/threads.h"
58 #include "gromacs/awh/awh.h"
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/compat/make_unique.h"
61 #include "gromacs/domdec/collect.h"
62 #include "gromacs/domdec/domdec.h"
63 #include "gromacs/domdec/domdec_network.h"
64 #include "gromacs/domdec/domdec_struct.h"
65 #include "gromacs/essentialdynamics/edsam.h"
66 #include "gromacs/ewald/pme.h"
67 #include "gromacs/ewald/pme-load-balancing.h"
68 #include "gromacs/fileio/trxio.h"
69 #include "gromacs/gmxlib/network.h"
70 #include "gromacs/gmxlib/nrnb.h"
71 #include "gromacs/gpu_utils/gpu_utils.h"
72 #include "gromacs/imd/imd.h"
73 #include "gromacs/listed-forces/manage-threading.h"
74 #include "gromacs/math/functions.h"
75 #include "gromacs/math/utilities.h"
76 #include "gromacs/math/vec.h"
77 #include "gromacs/math/vectypes.h"
78 #include "gromacs/mdlib/compute_io.h"
79 #include "gromacs/mdlib/constr.h"
80 #include "gromacs/mdlib/deform.h"
81 #include "gromacs/mdlib/ebin.h"
82 #include "gromacs/mdlib/expanded.h"
83 #include "gromacs/mdlib/force.h"
84 #include "gromacs/mdlib/force_flags.h"
85 #include "gromacs/mdlib/forcerec.h"
86 #include "gromacs/mdlib/md_support.h"
87 #include "gromacs/mdlib/mdatoms.h"
88 #include "gromacs/mdlib/mdebin.h"
89 #include "gromacs/mdlib/mdoutf.h"
90 #include "gromacs/mdlib/mdrun.h"
91 #include "gromacs/mdlib/mdsetup.h"
92 #include "gromacs/mdlib/membed.h"
93 #include "gromacs/mdlib/nb_verlet.h"
94 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
95 #include "gromacs/mdlib/ns.h"
96 #include "gromacs/mdlib/repl_ex.h"
97 #include "gromacs/mdlib/shellfc.h"
98 #include "gromacs/mdlib/sighandler.h"
99 #include "gromacs/mdlib/sim_util.h"
100 #include "gromacs/mdlib/simulationsignal.h"
101 #include "gromacs/mdlib/tgroup.h"
102 #include "gromacs/mdlib/trajectory_writing.h"
103 #include "gromacs/mdlib/update.h"
104 #include "gromacs/mdlib/vcm.h"
105 #include "gromacs/mdlib/vsite.h"
106 #include "gromacs/mdtypes/awh-history.h"
107 #include "gromacs/mdtypes/awh-params.h"
108 #include "gromacs/mdtypes/commrec.h"
109 #include "gromacs/mdtypes/df_history.h"
110 #include "gromacs/mdtypes/energyhistory.h"
111 #include "gromacs/mdtypes/fcdata.h"
112 #include "gromacs/mdtypes/forcerec.h"
113 #include "gromacs/mdtypes/group.h"
114 #include "gromacs/mdtypes/inputrec.h"
115 #include "gromacs/mdtypes/interaction_const.h"
116 #include "gromacs/mdtypes/md_enums.h"
117 #include "gromacs/mdtypes/mdatom.h"
118 #include "gromacs/mdtypes/observableshistory.h"
119 #include "gromacs/mdtypes/state.h"
120 #include "gromacs/pbcutil/mshift.h"
121 #include "gromacs/pbcutil/pbc.h"
122 #include "gromacs/pulling/pull.h"
123 #include "gromacs/swap/swapcoords.h"
124 #include "gromacs/timing/wallcycle.h"
125 #include "gromacs/timing/walltime_accounting.h"
126 #include "gromacs/topology/atoms.h"
127 #include "gromacs/topology/idef.h"
128 #include "gromacs/topology/mtop_util.h"
129 #include "gromacs/topology/topology.h"
130 #include "gromacs/trajectory/trajectoryframe.h"
131 #include "gromacs/utility/basedefinitions.h"
132 #include "gromacs/utility/cstringutil.h"
133 #include "gromacs/utility/fatalerror.h"
134 #include "gromacs/utility/logger.h"
135 #include "gromacs/utility/real.h"
136 #include "gromacs/utility/smalloc.h"
138 #include "integrator.h"
141 #include "corewrap.h"
144 using gmx::SimulationSignaller;
146 //! Resets all the counters.
147 static void reset_all_counters(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
149 gmx_int64_t *step_rel, t_inputrec *ir,
150 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
151 gmx_walltime_accounting_t walltime_accounting,
152 struct nonbonded_verlet_t *nbv,
153 struct gmx_pme_t *pme)
155 char sbuf[STEPSTRSIZE];
157 /* Reset all the counters related to performance over the run */
158 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
159 "step %s: resetting all time and cycle counters",
160 gmx_step_str(step, sbuf));
164 nbnxn_gpu_reset_timings(nbv);
167 if (pme_gpu_task_enabled(pme))
169 pme_gpu_reset_timings(pme);
172 if (use_GPU(nbv) || pme_gpu_task_enabled(pme))
177 wallcycle_stop(wcycle, ewcRUN);
178 wallcycle_reset_all(wcycle);
179 if (DOMAINDECOMP(cr))
181 reset_dd_statistics_counters(cr->dd);
184 ir->init_step += *step_rel;
185 ir->nsteps -= *step_rel;
187 wallcycle_start(wcycle, ewcRUN);
188 walltime_accounting_start(walltime_accounting);
189 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
192 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
194 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
195 * \param[in,out] globalState The global state container
196 * \param[in] constructVsites When true, vsite coordinates are constructed
197 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
198 * \param[in] idef Topology parameters, used for constructing vsites
199 * \param[in] timeStep Time step, used for constructing vsites
200 * \param[in] forceRec Force record, used for constructing vsites
201 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
202 * \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
204 static void prepareRerunState(const t_trxframe &rerunFrame,
205 t_state *globalState,
206 bool constructVsites,
207 const gmx_vsite_t *vsite,
210 const t_forcerec &forceRec,
212 gmx_bool *warnWhenNoV)
214 for (int i = 0; i < globalState->natoms; i++)
216 copy_rvec(rerunFrame.x[i], globalState->x[i]);
220 for (int i = 0; i < globalState->natoms; i++)
222 copy_rvec(rerunFrame.v[i], globalState->v[i]);
227 for (int i = 0; i < globalState->natoms; i++)
229 clear_rvec(globalState->v[i]);
233 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
234 " Ekin, temperature and pressure are incorrect,\n"
235 " the virial will be incorrect when constraints are present.\n"
237 *warnWhenNoV = FALSE;
240 copy_mat(rerunFrame.box, globalState->box);
244 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
248 /* Following is necessary because the graph may get out of sync
249 * with the coordinates if we only have every N'th coordinate set
251 mk_mshift(nullptr, graph, forceRec.ePBC, globalState->box, as_rvec_array(globalState->x.data()));
252 shift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
254 construct_vsites(vsite, as_rvec_array(globalState->x.data()), timeStep, as_rvec_array(globalState->v.data()),
255 idef.iparams, idef.il,
256 forceRec.ePBC, forceRec.bMolPBC, nullptr, globalState->box);
259 unshift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
264 void gmx::Integrator::do_md()
266 // TODO Historically, the EM and MD "integrators" used different
267 // names for the t_inputrec *parameter, but these must have the
268 // same name, now that it's a member of a struct. We use this ir
269 // alias to avoid a large ripple of nearly useless changes.
270 // t_inputrec is being replaced by IMdpOptionsProvider, so this
271 // will go away eventually.
272 t_inputrec *ir = inputrec;
273 gmx_mdoutf *outf = nullptr;
274 gmx_int64_t step, step_rel;
276 double t, t0, lam0[efptNR];
277 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
278 gmx_bool bNS, bNStList, bSimAnn, bStopCM,
279 bFirstStep, bInitStep, bLastStep = FALSE;
280 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
281 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
282 bForceUpdate = FALSE, bCPT;
283 gmx_bool bMasterState;
284 int force_flags, cglo_flags;
285 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
290 matrix parrinellorahmanMu, M;
292 gmx_repl_ex_t repl_ex = nullptr;
295 t_mdebin *mdebin = nullptr;
296 gmx_enerdata_t *enerd;
297 PaddedRVecVector f {};
298 gmx_global_stat_t gstat;
299 gmx_update_t *upd = nullptr;
300 t_graph *graph = nullptr;
301 gmx_groups_t *groups;
302 gmx_ekindata_t *ekind;
303 gmx_shellfc_t *shellfc;
304 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
305 gmx_bool bResetCountersHalfMaxH = FALSE;
306 gmx_bool bTemp, bPres, bTrotter;
308 rvec *cbuf = nullptr;
315 real saved_conserved_quantity = 0;
319 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
320 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
323 /* PME load balancing data for GPU kernels */
324 pme_load_balancing_t *pme_loadbal = nullptr;
325 gmx_bool bPMETune = FALSE;
326 gmx_bool bPMETunePrinting = FALSE;
329 gmx_bool bIMDstep = FALSE;
332 /* Temporary addition for FAHCORE checkpointing */
335 /* Domain decomposition could incorrectly miss a bonded
336 interaction, but checking for that requires a global
337 communication stage, which does not otherwise happen in DD
338 code. So we do that alongside the first global energy reduction
339 after a new DD is made. These variables handle whether the
340 check happens, and the result it returns. */
341 bool shouldCheckNumberOfBondedInteractions = false;
342 int totalNumberOfBondedInteractions = -1;
344 SimulationSignals signals;
345 // Most global communnication stages don't propagate mdrun
346 // signals, and will use this object to achieve that.
347 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
349 if (!mdrunOptions.writeConfout)
351 // This is on by default, and the main known use case for
352 // turning it off is for convenience in benchmarking, which is
353 // something that should not show up in the general user
355 GMX_LOG(mdlog.info).asParagraph().
356 appendText("The -noconfout functionality is deprecated, and may be removed in a future version.");
358 if (mdrunOptions.timingOptions.resetHalfway)
360 GMX_LOG(mdlog.info).asParagraph().
361 appendText("The -resethway functionality is deprecated, and may be removed in a future version.");
364 /* Signal to reset the counters half the simulation steps. */
365 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
367 /* Signal to reset the counters halfway the simulation time. */
368 bResetCountersHalfMaxH = (mdrunOptions.maximumHoursToRun > 0);
371 /* md-vv uses averaged full step velocities for T-control
372 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
373 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
374 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
376 const gmx_bool bRerunMD = mdrunOptions.rerun;
377 int nstglobalcomm = mdrunOptions.globalCommunicationInterval;
381 /* Since we don't know if the frames read are related in any way,
382 * rebuild the neighborlist at every step.
385 ir->nstcalcenergy = 1;
389 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
390 bGStatEveryStep = (nstglobalcomm == 1);
394 ir->nstxout_compressed = 0;
396 groups = &top_global->groups;
398 gmx_edsam *ed = nullptr;
399 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
401 /* Initialize essential dynamics sampling */
402 ed = init_edsam(opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm),
405 state_global, observablesHistory,
406 oenv, mdrunOptions.continuationOptions.appendFiles);
409 if (ir->eSwapCoords != eswapNO)
411 /* Initialize ion swapping code */
412 init_swapcoords(fplog, ir, opt2fn_master("-swap", nfile, fnm, cr),
414 state_global, observablesHistory,
415 cr, oenv, mdrunOptions);
419 init_md(fplog, cr, outputProvider, ir, oenv, mdrunOptions,
420 &t, &t0, state_global, lam0,
421 nrnb, top_global, &upd,
422 nfile, fnm, &outf, &mdebin,
423 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, wcycle);
425 clear_mat(total_vir);
427 /* Energy terms and groups */
429 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
432 /* Kinetic energy data */
434 init_ekindata(fplog, top_global, &(ir->opts), ekind);
435 /* Copy the cos acceleration to the groups struct */
436 ekind->cosacc.cos_accel = ir->cos_accel;
438 gstat = global_stat_init(ir);
440 /* Check for polarizable models and flexible constraints */
441 shellfc = init_shell_flexcon(fplog,
442 top_global, n_flexible_constraints(constr),
443 ir->nstcalcenergy, DOMAINDECOMP(cr));
445 if (shellfc && ir->nstcalcenergy != 1)
447 gmx_fatal(FARGS, "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is not supported in combinations with shell particles.\nPlease make a new tpr file.", ir->nstcalcenergy);
449 if (shellfc && DOMAINDECOMP(cr))
451 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
453 if (shellfc && ir->bDoAwh)
455 gmx_fatal(FARGS, "AWH biasing does not support shell particles.");
458 if (inputrecDeform(ir))
460 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
461 set_deform_reference_box(upd,
462 deform_init_init_step_tpx,
463 deform_init_box_tpx);
464 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
468 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
469 if ((io > 2000) && MASTER(cr))
472 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
477 // Local state only becomes valid now.
478 std::unique_ptr<t_state> stateInstance;
481 if (DOMAINDECOMP(cr))
483 top = dd_init_local_top(top_global);
485 stateInstance = gmx::compat::make_unique<t_state>();
486 state = stateInstance.get();
487 dd_init_local_state(cr->dd, state_global, state);
491 state_change_natoms(state_global, state_global->natoms);
492 /* We need to allocate one element extra, since we might use
493 * (unaligned) 4-wide SIMD loads to access rvec entries.
495 f.resize(gmx::paddedRVecVectorSize(state_global->natoms));
496 /* Copy the pointer to the global state */
497 state = state_global;
500 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
501 &graph, mdAtoms, vsite, shellfc);
503 update_realloc(upd, state->natoms);
506 /* Set up interactive MD (IMD) */
507 init_IMD(ir, cr, ms, top_global, fplog, ir->nstcalcenergy,
508 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
509 nfile, fnm, oenv, mdrunOptions);
511 if (DOMAINDECOMP(cr))
513 /* Distribute the charge groups over the nodes from the master node */
514 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
515 state_global, top_global, ir,
516 state, &f, mdAtoms, top, fr,
518 nrnb, nullptr, FALSE);
519 shouldCheckNumberOfBondedInteractions = true;
520 update_realloc(upd, state->natoms);
523 auto mdatoms = mdAtoms->mdatoms();
525 // NOTE: The global state is no longer used at this point.
526 // But state_global is still used as temporary storage space for writing
527 // the global state to file and potentially for replica exchange.
528 // (Global topology should persist.)
530 update_mdatoms(mdatoms, state->lambda[efptMASS]);
532 const ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
533 bool startingFromCheckpoint = continuationOptions.startedFromCheckpoint;
537 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
542 if (startingFromCheckpoint)
544 /* Update mdebin with energy history if appending to output files */
545 if (continuationOptions.appendFiles)
547 restore_energyhistory_from_state(mdebin, observablesHistory->energyHistory.get());
549 else if (observablesHistory->energyHistory.get() != nullptr)
551 /* We might have read an energy history from checkpoint.
552 * As we are not appending, we want to restart the statistics.
553 * Free the allocated memory and reset the counts.
555 observablesHistory->energyHistory = {};
558 if (observablesHistory->energyHistory.get() == nullptr)
560 observablesHistory->energyHistory = std::unique_ptr<energyhistory_t>(new energyhistory_t {});
562 /* Set the initial energy history in state by updating once */
563 update_energyhistory(observablesHistory->energyHistory.get(), mdebin);
566 /* Initialize constraints */
567 if (constr && !DOMAINDECOMP(cr))
569 set_constraints(constr, top, ir, mdatoms, cr);
572 /* Initialize AWH and restore state from history in checkpoint if needed. */
575 ir->awh = new gmx::Awh(fplog, *ir, cr, ms, *ir->awhParams, opt2fn("-awh", nfile, fnm), ir->pull_work);
577 if (startingFromCheckpoint)
579 /* Restore the AWH history read from checkpoint */
580 ir->awh->restoreStateFromHistory(MASTER(cr) ? state_global->awhHistory.get() : nullptr);
584 /* Initialize the AWH history here */
585 state_global->awhHistory = ir->awh->initHistoryFromState();
589 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
590 if (useReplicaExchange && MASTER(cr))
592 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir,
595 /* PME tuning is only supported in the Verlet scheme, with PME for
596 * Coulomb. It is not supported with only LJ PME, or for
598 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !bRerunMD &&
599 !mdrunOptions.reproducible && ir->cutoff_scheme != ecutsGROUP);
602 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
603 fr->ic, fr->nbv->listParams.get(), fr->pmedata, use_GPU(fr->nbv),
607 if (!ir->bContinuation && !bRerunMD)
609 if (state->flags & (1 << estV))
611 /* Set the velocities of vsites, shells and frozen atoms to zero */
612 for (i = 0; i < mdatoms->homenr; i++)
614 if (mdatoms->ptype[i] == eptVSite ||
615 mdatoms->ptype[i] == eptShell)
617 clear_rvec(state->v[i]);
619 else if (mdatoms->cFREEZE)
621 for (m = 0; m < DIM; m++)
623 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
634 /* Constrain the initial coordinates and velocities */
635 do_constrain_first(fplog, constr, ir, mdatoms, state,
636 cr, ms, nrnb, fr, top);
640 /* Construct the virtual sites for the initial configuration */
641 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, nullptr,
642 top->idef.iparams, top->idef.il,
643 fr->ePBC, fr->bMolPBC, cr, state->box);
647 if (ir->efep != efepNO)
649 /* Set free energy calculation frequency as the greatest common
650 * denominator of nstdhdl and repl_ex_nst. */
651 nstfep = ir->fepvals->nstdhdl;
654 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
656 if (useReplicaExchange)
658 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
662 /* Be REALLY careful about what flags you set here. You CANNOT assume
663 * this is the first step, since we might be restarting from a checkpoint,
664 * and in that case we should not do any modifications to the state.
666 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
668 if (continuationOptions.haveReadEkin)
670 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
673 cglo_flags = (CGLO_INITIALIZATION | CGLO_TEMPERATURE | CGLO_GSTAT
674 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
675 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
676 | (continuationOptions.haveReadEkin ? CGLO_READEKIN : 0));
678 bSumEkinhOld = FALSE;
679 /* To minimize communication, compute_globals computes the COM velocity
680 * and the kinetic energy for the velocities without COM motion removed.
681 * Thus to get the kinetic energy without the COM contribution, we need
682 * to call compute_globals twice.
684 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
686 int cglo_flags_iteration = cglo_flags;
687 if (bStopCM && cgloIteration == 0)
689 cglo_flags_iteration |= CGLO_STOPCM;
690 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
692 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
693 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
694 constr, &nullSignaller, state->box,
695 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags_iteration
696 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
698 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
699 top_global, top, state,
700 &shouldCheckNumberOfBondedInteractions);
701 if (ir->eI == eiVVAK)
703 /* a second call to get the half step temperature initialized as well */
704 /* we do the same call as above, but turn the pressure off -- internally to
705 compute_globals, this is recognized as a velocity verlet half-step
706 kinetic energy calculation. This minimized excess variables, but
707 perhaps loses some logic?*/
709 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
710 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
711 constr, &nullSignaller, state->box,
712 nullptr, &bSumEkinhOld,
713 cglo_flags & ~CGLO_PRESSURE);
716 /* Calculate the initial half step temperature, and save the ekinh_old */
717 if (!continuationOptions.startedFromCheckpoint)
719 for (i = 0; (i < ir->opts.ngtc); i++)
721 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
725 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
726 temperature control */
727 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
731 if (!ir->bContinuation)
733 if (constr && ir->eConstrAlg == econtLINCS)
736 "RMS relative constraint deviation after constraining: %.2e\n",
737 constr_rmsd(constr));
739 if (EI_STATE_VELOCITY(ir->eI))
741 real temp = enerd->term[F_TEMP];
744 /* Result of Ekin averaged over velocities of -half
745 * and +half step, while we only have -half step here.
749 fprintf(fplog, "Initial temperature: %g K\n", temp);
755 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
756 " input trajectory '%s'\n\n",
757 *(top_global->name), opt2fn("-rerun", nfile, fnm));
758 if (mdrunOptions.verbose)
760 fprintf(stderr, "Calculated time to finish depends on nsteps from "
761 "run input file,\nwhich may not correspond to the time "
762 "needed to process input trajectory.\n\n");
768 fprintf(stderr, "starting mdrun '%s'\n",
769 *(top_global->name));
772 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
776 sprintf(tbuf, "%s", "infinite");
778 if (ir->init_step > 0)
780 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
781 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
782 gmx_step_str(ir->init_step, sbuf2),
783 ir->init_step*ir->delta_t);
787 fprintf(stderr, "%s steps, %s ps.\n",
788 gmx_step_str(ir->nsteps, sbuf), tbuf);
791 fprintf(fplog, "\n");
794 walltime_accounting_start(walltime_accounting);
795 wallcycle_start(wcycle, ewcRUN);
796 print_start(fplog, cr, walltime_accounting, "mdrun");
798 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
800 chkpt_ret = fcCheckPointParallel( cr->nodeid,
804 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
808 /***********************************************************
812 ************************************************************/
814 /* if rerunMD then read coordinates and velocities from input trajectory */
817 if (getenv("GMX_FORCE_UPDATE"))
825 bLastStep = !read_first_frame(oenv, &status,
826 opt2fn("-rerun", nfile, fnm),
827 &rerun_fr, TRX_NEED_X | TRX_READ_V);
828 if (rerun_fr.natoms != top_global->natoms)
831 "Number of atoms in trajectory (%d) does not match the "
832 "run input file (%d)\n",
833 rerun_fr.natoms, top_global->natoms);
835 if (ir->ePBC != epbcNONE)
839 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f does not contain a box, while pbc is used", rerun_fr.step, rerun_fr.time);
841 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
843 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
850 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
853 if (ir->ePBC != epbcNONE)
855 /* Set the shift vectors.
856 * Necessary here when have a static box different from the tpr box.
858 calc_shifts(rerun_fr.box, fr->shift_vec);
862 /* Loop over MD steps or if rerunMD to end of input trajectory,
863 * or, if max_hours>0, until max_hours is reached.
865 real max_hours = mdrunOptions.maximumHoursToRun;
867 /* Skip the first Nose-Hoover integration when we get the state from tpx */
868 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
869 bSumEkinhOld = FALSE;
871 bNeedRepartition = FALSE;
873 bool simulationsShareState = false;
874 int nstSignalComm = nstglobalcomm;
876 // TODO This implementation of ensemble orientation restraints is nasty because
877 // a user can't just do multi-sim with single-sim orientation restraints.
878 bool usingEnsembleRestraints = (fcd->disres.nsystems > 1) || (ms && fcd->orires.nr);
879 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && ms);
881 // Replica exchange, ensemble restraints and AWH need all
882 // simulations to remain synchronized, so they need
883 // checkpoints and stop conditions to act on the same step, so
884 // the propagation of such signals must take place between
885 // simulations, not just within simulations.
886 // TODO: Make algorithm initializers set these flags.
887 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
888 bool resetCountersIsLocal = true;
889 signals[eglsCHKPT] = SimulationSignal(!simulationsShareState);
890 signals[eglsSTOPCOND] = SimulationSignal(!simulationsShareState);
891 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
893 if (simulationsShareState)
895 // Inter-simulation signal communication does not need to happen
896 // often, so we use a minimum of 200 steps to reduce overhead.
897 const int c_minimumInterSimulationSignallingInterval = 200;
898 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1)/nstglobalcomm)*nstglobalcomm;
902 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion = (DOMAINDECOMP(cr) ? DdOpenBalanceRegionBeforeForceComputation::yes : DdOpenBalanceRegionBeforeForceComputation::no);
903 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion = (DOMAINDECOMP(cr) ? DdCloseBalanceRegionAfterForceComputation::yes : DdCloseBalanceRegionAfterForceComputation::no);
905 step = ir->init_step;
908 // TODO extract this to new multi-simulation module
909 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
911 if (!multisim_int_all_are_equal(ms, ir->nsteps))
913 GMX_LOG(mdlog.warning).appendText(
914 "Note: The number of steps is not consistent across multi simulations,\n"
915 "but we are proceeding anyway!");
917 if (!multisim_int_all_are_equal(ms, ir->init_step))
919 GMX_LOG(mdlog.warning).appendText(
920 "Note: The initial step is not consistent across multi simulations,\n"
921 "but we are proceeding anyway!");
925 /* and stop now if we should */
926 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
930 /* Determine if this is a neighbor search step */
931 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
933 if (bPMETune && bNStList)
935 /* PME grid + cut-off optimization with GPUs or PME nodes */
936 pme_loadbal_do(pme_loadbal, cr,
937 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
945 wallcycle_start(wcycle, ewcSTEP);
951 step = rerun_fr.step;
952 step_rel = step - ir->init_step;
965 bLastStep = (step_rel == ir->nsteps);
966 t = t0 + step*ir->delta_t;
969 // TODO Refactor this, so that nstfep does not need a default value of zero
970 if (ir->efep != efepNO || ir->bSimTemp)
972 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
973 requiring different logic. */
978 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
983 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
985 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
986 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
987 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
988 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
991 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep &&
992 do_per_step(step, replExParams.exchangeInterval));
996 update_annealing_target_temp(ir, t, upd);
999 if (bRerunMD && MASTER(cr))
1001 const bool constructVsites = (vsite && mdrunOptions.rerunConstructVsites);
1002 if (constructVsites && DOMAINDECOMP(cr))
1004 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
1006 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top->idef, ir->delta_t, *fr, graph, &bRerunWarnNoV);
1009 /* Stop Center of Mass motion */
1010 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
1014 /* for rerun MD always do Neighbour Searching */
1015 bNS = (bFirstStep || ir->nstlist != 0);
1019 /* Determine whether or not to do Neighbour Searching */
1020 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
1023 /* < 0 means stop at next step, > 0 means stop at next NS step */
1024 if ( (signals[eglsSTOPCOND].set < 0) ||
1025 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
1030 /* do_log triggers energy and virial calculation. Because this leads
1031 * to different code paths, forces can be different. Thus for exact
1032 * continuation we should avoid extra log output.
1033 * Note that the || bLastStep can result in non-exact continuation
1034 * beyond the last step. But we don't consider that to be an issue.
1036 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
1037 do_verbose = mdrunOptions.verbose &&
1038 (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep || bRerunMD);
1040 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
1044 bMasterState = TRUE;
1048 bMasterState = FALSE;
1049 /* Correct the new box if it is too skewed */
1050 if (inputrecDynamicBox(ir))
1052 if (correct_box(fplog, step, state->box, graph))
1054 bMasterState = TRUE;
1057 if (DOMAINDECOMP(cr) && bMasterState)
1059 dd_collect_state(cr->dd, state, state_global);
1063 if (DOMAINDECOMP(cr))
1065 /* Repartition the domain decomposition */
1066 dd_partition_system(fplog, step, cr,
1067 bMasterState, nstglobalcomm,
1068 state_global, top_global, ir,
1069 state, &f, mdAtoms, top, fr,
1072 do_verbose && !bPMETunePrinting);
1073 shouldCheckNumberOfBondedInteractions = true;
1074 update_realloc(upd, state->natoms);
1078 if (MASTER(cr) && do_log)
1080 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
1083 if (ir->efep != efepNO)
1085 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1088 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1091 /* We need the kinetic energy at minus the half step for determining
1092 * the full step kinetic energy and possibly for T-coupling.*/
1093 /* This may not be quite working correctly yet . . . . */
1094 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1095 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1096 constr, &nullSignaller, state->box,
1097 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1098 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1099 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1100 top_global, top, state,
1101 &shouldCheckNumberOfBondedInteractions);
1103 clear_mat(force_vir);
1105 /* We write a checkpoint at this MD step when:
1106 * either at an NS step when we signalled through gs,
1107 * or at the last step (but not when we do not want confout),
1108 * but never at the first step or with rerun.
1110 bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1111 (bLastStep && mdrunOptions.writeConfout)) &&
1112 step > ir->init_step && !bRerunMD);
1115 signals[eglsCHKPT].set = 0;
1118 /* Determine the energy and pressure:
1119 * at nstcalcenergy steps and at energy output steps (set below).
1121 if (EI_VV(ir->eI) && (!bInitStep))
1123 /* for vv, the first half of the integration actually corresponds
1124 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1125 but the virial needs to be calculated on both the current step and the 'next' step. Future
1126 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1128 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1129 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1130 bCalcVir = bCalcEnerStep ||
1131 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1135 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1136 bCalcVir = bCalcEnerStep ||
1137 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1139 bCalcEner = bCalcEnerStep;
1141 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1143 if (do_ene || do_log || bDoReplEx)
1149 /* Do we need global communication ? */
1150 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1151 do_per_step(step, nstglobalcomm) ||
1152 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1154 force_flags = (GMX_FORCE_STATECHANGED |
1155 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1156 GMX_FORCE_ALLFORCES |
1157 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1158 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1159 (bDoFEP ? GMX_FORCE_DHDL : 0)
1164 /* Now is the time to relax the shells */
1165 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, step,
1166 ir, bNS, force_flags, top,
1168 state, f, force_vir, mdatoms,
1169 nrnb, wcycle, graph, groups,
1170 shellfc, fr, t, mu_tot,
1172 ddOpenBalanceRegion, ddCloseBalanceRegion);
1176 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
1177 (or the AWH update will be performed twice for one step when continuing). It would be best to
1178 call this update function from do_md_trajectory_writing but that would occur after do_force.
1179 One would have to divide the update_awh function into one function applying the AWH force
1180 and one doing the AWH bias update. The update AWH bias function could then be called after
1181 do_md_trajectory_writing (then containing update_awh_history).
1182 The checkpointing will in the future probably moved to the start of the md loop which will
1183 rid of this issue. */
1184 if (ir->bDoAwh && bCPT && MASTER(cr))
1186 ir->awh->updateHistory(state_global->awhHistory.get());
1189 /* The coordinates (x) are shifted (to get whole molecules)
1191 * This is parallellized as well, and does communication too.
1192 * Check comments in sim_util.c
1194 do_force(fplog, cr, ms, ir, step, nrnb, wcycle, top, groups,
1195 state->box, state->x, &state->hist,
1196 f, force_vir, mdatoms, enerd, fcd,
1197 state->lambda, graph,
1198 fr, vsite, mu_tot, t, ed,
1199 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1200 ddOpenBalanceRegion, ddCloseBalanceRegion);
1203 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1204 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1206 rvec *vbuf = nullptr;
1208 wallcycle_start(wcycle, ewcUPDATE);
1209 if (ir->eI == eiVV && bInitStep)
1211 /* if using velocity verlet with full time step Ekin,
1212 * take the first half step only to compute the
1213 * virial for the first step. From there,
1214 * revert back to the initial coordinates
1215 * so that the input is actually the initial step.
1217 snew(vbuf, state->natoms);
1218 copy_rvecn(as_rvec_array(state->v.data()), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1222 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1223 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1226 update_coords(step, ir, mdatoms, state, f, fcd,
1227 ekind, M, upd, etrtVELOCITY1,
1230 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1232 wallcycle_stop(wcycle, ewcUPDATE);
1233 constrain_velocities(step, nullptr, ir, mdatoms,
1235 &top->idef, shake_vir,
1236 cr, ms, nrnb, wcycle, constr,
1237 bCalcVir, do_log, do_ene);
1238 wallcycle_start(wcycle, ewcUPDATE);
1242 /* Need to unshift here if a do_force has been
1243 called in the previous step */
1244 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1246 /* if VV, compute the pressure and constraints */
1247 /* For VV2, we strictly only need this if using pressure
1248 * control, but we really would like to have accurate pressures
1250 * Think about ways around this in the future?
1251 * For now, keep this choice in comments.
1253 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1254 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1256 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1257 if (bCalcEner && ir->eI == eiVVAK)
1259 bSumEkinhOld = TRUE;
1261 /* for vv, the first half of the integration actually corresponds to the previous step.
1262 So we need information from the last step in the first half of the integration */
1263 if (bGStat || do_per_step(step-1, nstglobalcomm))
1265 wallcycle_stop(wcycle, ewcUPDATE);
1266 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1267 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1268 constr, &nullSignaller, state->box,
1269 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1270 (bGStat ? CGLO_GSTAT : 0)
1272 | (bTemp ? CGLO_TEMPERATURE : 0)
1273 | (bPres ? CGLO_PRESSURE : 0)
1274 | (bPres ? CGLO_CONSTRAINT : 0)
1275 | (bStopCM ? CGLO_STOPCM : 0)
1276 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1279 /* explanation of above:
1280 a) We compute Ekin at the full time step
1281 if 1) we are using the AveVel Ekin, and it's not the
1282 initial step, or 2) if we are using AveEkin, but need the full
1283 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1284 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1285 EkinAveVel because it's needed for the pressure */
1286 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1287 top_global, top, state,
1288 &shouldCheckNumberOfBondedInteractions);
1289 wallcycle_start(wcycle, ewcUPDATE);
1291 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1296 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1297 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1299 /* TODO This is only needed when we're about to write
1300 * a checkpoint, because we use it after the restart
1301 * (in a kludge?). But what should we be doing if
1302 * startingFromCheckpoint or bInitStep are true? */
1303 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1305 copy_mat(shake_vir, state->svir_prev);
1306 copy_mat(force_vir, state->fvir_prev);
1308 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1310 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1311 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1312 enerd->term[F_EKIN] = trace(ekind->ekin);
1315 else if (bExchanged)
1317 wallcycle_stop(wcycle, ewcUPDATE);
1318 /* We need the kinetic energy at minus the half step for determining
1319 * the full step kinetic energy and possibly for T-coupling.*/
1320 /* This may not be quite working correctly yet . . . . */
1321 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1322 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1323 constr, &nullSignaller, state->box,
1324 nullptr, &bSumEkinhOld,
1325 CGLO_GSTAT | CGLO_TEMPERATURE);
1326 wallcycle_start(wcycle, ewcUPDATE);
1329 /* if it's the initial step, we performed this first step just to get the constraint virial */
1330 if (ir->eI == eiVV && bInitStep)
1332 copy_rvecn(vbuf, as_rvec_array(state->v.data()), 0, state->natoms);
1335 wallcycle_stop(wcycle, ewcUPDATE);
1338 /* compute the conserved quantity */
1341 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1344 last_ekin = enerd->term[F_EKIN];
1346 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1348 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1350 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1351 if (ir->efep != efepNO && !bRerunMD)
1353 sum_dhdl(enerd, state->lambda, ir->fepvals);
1357 /* ######## END FIRST UPDATE STEP ############## */
1358 /* ######## If doing VV, we now have v(dt) ###### */
1361 /* perform extended ensemble sampling in lambda - we don't
1362 actually move to the new state before outputting
1363 statistics, but if performing simulated tempering, we
1364 do update the velocities and the tau_t. */
1366 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, as_rvec_array(state->v.data()), mdatoms);
1367 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1370 copy_df_history(state_global->dfhist, state->dfhist);
1374 /* Now we have the energies and forces corresponding to the
1375 * coordinates at time t. We must output all of this before
1378 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1379 ir, state, state_global, observablesHistory,
1381 outf, mdebin, ekind, f,
1383 bCPT, bRerunMD, bLastStep,
1384 mdrunOptions.writeConfout,
1386 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1387 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, as_rvec_array(state->x.data()), ir, t, wcycle);
1389 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1390 if (startingFromCheckpoint && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1392 copy_mat(state->svir_prev, shake_vir);
1393 copy_mat(state->fvir_prev, force_vir);
1396 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1398 /* Check whether everything is still allright */
1399 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1405 int nsteps_stop = -1;
1407 /* this just makes signals[].sig compatible with the hack
1408 of sending signals around by MPI_Reduce together with
1410 if ((gmx_get_stop_condition() == gmx_stop_cond_next_ns) ||
1411 (mdrunOptions.reproducible &&
1412 gmx_get_stop_condition() == gmx_stop_cond_next))
1414 /* We need at least two global communication steps to pass
1415 * around the signal. We stop at a pair-list creation step
1416 * to allow for exact continuation, when possible.
1418 signals[eglsSTOPCOND].sig = 1;
1419 nsteps_stop = std::max(ir->nstlist, 2*nstSignalComm);
1421 else if (gmx_get_stop_condition() == gmx_stop_cond_next)
1423 /* Stop directly after the next global communication step.
1424 * This breaks exact continuation.
1426 signals[eglsSTOPCOND].sig = -1;
1427 nsteps_stop = nstSignalComm + 1;
1432 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1433 gmx_get_signal_name(), nsteps_stop);
1437 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1438 gmx_get_signal_name(), nsteps_stop);
1440 handled_stop_condition = (int)gmx_get_stop_condition();
1442 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1443 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1444 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1446 /* Signal to terminate the run */
1447 signals[eglsSTOPCOND].sig = 1;
1450 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1452 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1455 if (bResetCountersHalfMaxH && MASTER(cr) &&
1456 elapsed_time > max_hours*60.0*60.0*0.495)
1458 /* Set flag that will communicate the signal to all ranks in the simulation */
1459 signals[eglsRESETCOUNTERS].sig = 1;
1462 /* In parallel we only have to check for checkpointing in steps
1463 * where we do global communication,
1464 * otherwise the other nodes don't know.
1466 const real cpt_period = mdrunOptions.checkpointOptions.period;
1467 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1470 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1471 signals[eglsCHKPT].set == 0)
1473 signals[eglsCHKPT].sig = 1;
1476 /* ######### START SECOND UPDATE STEP ################# */
1478 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1481 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1483 gmx_bool bIfRandomize;
1484 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1485 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1486 if (constr && bIfRandomize)
1488 constrain_velocities(step, nullptr, ir, mdatoms,
1490 &top->idef, tmp_vir,
1491 cr, ms, nrnb, wcycle, constr,
1492 bCalcVir, do_log, do_ene);
1495 /* Box is changed in update() when we do pressure coupling,
1496 * but we should still use the old box for energy corrections and when
1497 * writing it to the energy file, so it matches the trajectory files for
1498 * the same timestep above. Make a copy in a separate array.
1500 copy_mat(state->box, lastbox);
1504 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1506 wallcycle_start(wcycle, ewcUPDATE);
1507 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1510 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1511 /* We can only do Berendsen coupling after we have summed
1512 * the kinetic energy or virial. Since the happens
1513 * in global_state after update, we should only do it at
1514 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1519 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1520 update_pcouple_before_coordinates(fplog, step, ir, state,
1521 parrinellorahmanMu, M,
1527 /* velocity half-step update */
1528 update_coords(step, ir, mdatoms, state, f, fcd,
1529 ekind, M, upd, etrtVELOCITY2,
1533 /* Above, initialize just copies ekinh into ekin,
1534 * it doesn't copy position (for VV),
1535 * and entire integrator for MD.
1538 if (ir->eI == eiVVAK)
1540 /* We probably only need md->homenr, not state->natoms */
1541 if (state->natoms > cbuf_nalloc)
1543 cbuf_nalloc = state->natoms;
1544 srenew(cbuf, cbuf_nalloc);
1546 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1549 update_coords(step, ir, mdatoms, state, f, fcd,
1550 ekind, M, upd, etrtPOSITION, cr, constr);
1551 wallcycle_stop(wcycle, ewcUPDATE);
1553 constrain_coordinates(step, &dvdl_constr, ir, mdatoms, state,
1555 &top->idef, shake_vir,
1556 cr, ms, nrnb, wcycle, upd, constr,
1557 bCalcVir, do_log, do_ene);
1558 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state,
1559 fr->bMolPBC, &top->idef, cr, ms,
1560 nrnb, wcycle, upd, constr, do_log, do_ene);
1561 finish_update(ir, mdatoms,
1563 nrnb, wcycle, upd, constr);
1565 if (ir->eI == eiVVAK)
1567 /* erase F_EKIN and F_TEMP here? */
1568 /* just compute the kinetic energy at the half step to perform a trotter step */
1569 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1570 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1571 constr, &nullSignaller, lastbox,
1572 nullptr, &bSumEkinhOld,
1573 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1575 wallcycle_start(wcycle, ewcUPDATE);
1576 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1577 /* now we know the scaling, we can compute the positions again again */
1578 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1580 update_coords(step, ir, mdatoms, state, f, fcd,
1581 ekind, M, upd, etrtPOSITION, cr, constr);
1582 wallcycle_stop(wcycle, ewcUPDATE);
1584 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1585 /* are the small terms in the shake_vir here due
1586 * to numerical errors, or are they important
1587 * physically? I'm thinking they are just errors, but not completely sure.
1588 * For now, will call without actually constraining, constr=NULL*/
1589 finish_update(ir, mdatoms,
1591 nrnb, wcycle, upd, nullptr);
1595 /* this factor or 2 correction is necessary
1596 because half of the constraint force is removed
1597 in the vv step, so we have to double it. See
1598 the Redmine issue #1255. It is not yet clear
1599 if the factor of 2 is exact, or just a very
1600 good approximation, and this will be
1601 investigated. The next step is to see if this
1602 can be done adding a dhdl contribution from the
1603 rattle step, but this is somewhat more
1604 complicated with the current code. Will be
1605 investigated, hopefully for 4.6.3. However,
1606 this current solution is much better than
1607 having it completely wrong.
1609 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1613 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1618 /* Need to unshift here */
1619 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1622 if (vsite != nullptr)
1624 wallcycle_start(wcycle, ewcVSITECONSTR);
1625 if (graph != nullptr)
1627 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1629 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
1630 top->idef.iparams, top->idef.il,
1631 fr->ePBC, fr->bMolPBC, cr, state->box);
1633 if (graph != nullptr)
1635 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1637 wallcycle_stop(wcycle, ewcVSITECONSTR);
1640 /* ############## IF NOT VV, Calculate globals HERE ############ */
1641 /* With Leap-Frog we can skip compute_globals at
1642 * non-communication steps, but we need to calculate
1643 * the kinetic energy one step before communication.
1646 // Organize to do inter-simulation signalling on steps if
1647 // and when algorithms require it.
1648 bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1650 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1652 // Since we're already communicating at this step, we
1653 // can propagate intra-simulation signals. Note that
1654 // check_nstglobalcomm has the responsibility for
1655 // choosing the value of nstglobalcomm that is one way
1656 // bGStat becomes true, so we can't get into a
1657 // situation where e.g. checkpointing can't be
1659 bool doIntraSimSignal = true;
1660 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1662 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1663 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1666 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1667 (bGStat ? CGLO_GSTAT : 0)
1668 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1669 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1670 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1671 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1673 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1675 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1676 top_global, top, state,
1677 &shouldCheckNumberOfBondedInteractions);
1681 /* ############# END CALC EKIN AND PRESSURE ################# */
1683 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1684 the virial that should probably be addressed eventually. state->veta has better properies,
1685 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1686 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1688 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1690 /* Sum up the foreign energy and dhdl terms for md and sd.
1691 Currently done every step so that dhdl is correct in the .edr */
1692 sum_dhdl(enerd, state->lambda, ir->fepvals);
1695 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1696 pres, force_vir, shake_vir,
1700 /* ################# END UPDATE STEP 2 ################# */
1701 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1703 /* The coordinates (x) were unshifted in update */
1706 /* We will not sum ekinh_old,
1707 * so signal that we still have to do it.
1709 bSumEkinhOld = TRUE;
1714 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1716 /* use the directly determined last velocity, not actually the averaged half steps */
1717 if (bTrotter && ir->eI == eiVV)
1719 enerd->term[F_EKIN] = last_ekin;
1721 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1723 if (integratorHasConservedEnergyQuantity(ir))
1727 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1731 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1734 /* ######### END PREPARING EDR OUTPUT ########### */
1740 if (fplog && do_log && bDoExpanded)
1742 /* only needed if doing expanded ensemble */
1743 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1744 state_global->dfhist, state->fep_state, ir->nstlog, step);
1748 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1749 t, mdatoms->tmass, enerd, state,
1750 ir->fepvals, ir->expandedvals, lastbox,
1751 shake_vir, force_vir, total_vir, pres,
1752 ekind, mu_tot, constr);
1756 upd_mdebin_step(mdebin);
1759 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1760 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1762 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1764 eprNORMAL, mdebin, fcd, groups, &(ir->opts), ir->awh);
1768 pull_print_output(ir->pull_work, step, t);
1771 if (do_per_step(step, ir->nstlog))
1773 if (fflush(fplog) != 0)
1775 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1781 /* Have to do this part _after_ outputting the logfile and the edr file */
1782 /* Gets written into the state at the beginning of next loop*/
1783 state->fep_state = lamnew;
1785 /* Print the remaining wall clock time for the run */
1786 if (isMasterSimMasterRank(ms, cr) &&
1787 (do_verbose || gmx_got_usr_signal()) &&
1792 fprintf(stderr, "\n");
1794 print_time(stderr, walltime_accounting, step, ir, cr);
1797 /* Ion/water position swapping.
1798 * Not done in last step since trajectory writing happens before this call
1799 * in the MD loop and exchanges would be lost anyway. */
1800 bNeedRepartition = FALSE;
1801 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1802 do_per_step(step, ir->swap->nstswap))
1804 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1805 bRerunMD ? rerun_fr.x : as_rvec_array(state->x.data()),
1806 bRerunMD ? rerun_fr.box : state->box,
1807 MASTER(cr) && mdrunOptions.verbose,
1810 if (bNeedRepartition && DOMAINDECOMP(cr))
1812 dd_collect_state(cr->dd, state, state_global);
1816 /* Replica exchange */
1820 bExchanged = replica_exchange(fplog, cr, ms, repl_ex,
1821 state_global, enerd,
1825 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1827 dd_partition_system(fplog, step, cr, TRUE, 1,
1828 state_global, top_global, ir,
1829 state, &f, mdAtoms, top, fr,
1831 nrnb, wcycle, FALSE);
1832 shouldCheckNumberOfBondedInteractions = true;
1833 update_realloc(upd, state->natoms);
1838 startingFromCheckpoint = false;
1840 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1841 /* With all integrators, except VV, we need to retain the pressure
1842 * at the current step for coupling at the next step.
1844 if ((state->flags & (1<<estPRES_PREV)) &&
1846 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1848 /* Store the pressure in t_state for pressure coupling
1849 * at the next MD step.
1851 copy_mat(pres, state->pres_prev);
1854 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1856 if ( (membed != nullptr) && (!bLastStep) )
1858 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1865 /* read next frame from input trajectory */
1866 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1871 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1875 cycles = wallcycle_stop(wcycle, ewcSTEP);
1876 if (DOMAINDECOMP(cr) && wcycle)
1878 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1881 if (!bRerunMD || !rerun_fr.bStep)
1883 /* increase the MD step number */
1888 /* TODO make a counter-reset module */
1889 /* If it is time to reset counters, set a flag that remains
1890 true until counters actually get reset */
1891 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1892 signals[eglsRESETCOUNTERS].set != 0)
1894 if (pme_loadbal_is_active(pme_loadbal))
1896 /* Do not permit counter reset while PME load
1897 * balancing is active. The only purpose for resetting
1898 * counters is to measure reliable performance data,
1899 * and that can't be done before balancing
1902 * TODO consider fixing this by delaying the reset
1903 * until after load balancing completes,
1904 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1905 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1906 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1907 "resetting counters later in the run, e.g. with gmx "
1908 "mdrun -resetstep.", step);
1910 reset_all_counters(fplog, mdlog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1911 use_GPU(fr->nbv) ? fr->nbv : nullptr, fr->pmedata);
1912 wcycle_set_reset_counters(wcycle, -1);
1913 if (!thisRankHasDuty(cr, DUTY_PME))
1915 /* Tell our PME node to reset its counters */
1916 gmx_pme_send_resetcounters(cr, step);
1918 /* Correct max_hours for the elapsed time */
1919 max_hours -= elapsed_time/(60.0*60.0);
1920 /* If mdrun -maxh -resethway was active, it can only trigger once */
1921 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1922 /* Reset can only happen once, so clear the triggering flag. */
1923 signals[eglsRESETCOUNTERS].set = 0;
1926 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1927 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1930 /* End of main MD loop */
1932 /* Closing TNG files can include compressing data. Therefore it is good to do that
1933 * before stopping the time measurements. */
1934 mdoutf_tng_close(outf);
1936 /* Stop measuring walltime */
1937 walltime_accounting_end(walltime_accounting);
1939 if (bRerunMD && MASTER(cr))
1944 if (!thisRankHasDuty(cr, DUTY_PME))
1946 /* Tell the PME only node to finish */
1947 gmx_pme_send_finish(cr);
1952 if (ir->nstcalcenergy > 0 && !bRerunMD)
1954 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1955 eprAVER, mdebin, fcd, groups, &(ir->opts), ir->awh);
1958 done_mdebin(mdebin);
1963 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
1966 done_shellfc(fplog, shellfc, step_rel);
1968 if (useReplicaExchange && MASTER(cr))
1970 print_replica_exchange_statistics(fplog, repl_ex);
1978 // Clean up swapcoords
1979 if (ir->eSwapCoords != eswapNO)
1981 finish_swapcoords(ir->swap);
1984 /* Do essential dynamics cleanup if needed. Close .edo file */
1987 /* IMD cleanup, if bIMD is TRUE. */
1988 IMD_finalize(ir->bIMD, ir->imd);
1990 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1991 if (step_rel >= wcycle_get_reset_counters(wcycle) &&
1992 signals[eglsRESETCOUNTERS].set == 0 &&
1993 !bResetCountersHalfMaxH)
1995 walltime_accounting_set_valid_finish(walltime_accounting);
1998 destroy_enerdata(enerd);
2000 mdAlgorithmsTearDownAtomData(fr->bonded_threading, top);
2001 fr->bonded_threading = nullptr;