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,2019, 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
54 #include "gromacs/awh/awh.h"
55 #include "gromacs/commandline/filenm.h"
56 #include "gromacs/compat/make_unique.h"
57 #include "gromacs/domdec/collect.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_network.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/domdec/partition.h"
62 #include "gromacs/essentialdynamics/edsam.h"
63 #include "gromacs/ewald/pme.h"
64 #include "gromacs/ewald/pme-load-balancing.h"
65 #include "gromacs/fileio/trxio.h"
66 #include "gromacs/gmxlib/network.h"
67 #include "gromacs/gmxlib/nrnb.h"
68 #include "gromacs/gpu_utils/gpu_utils.h"
69 #include "gromacs/imd/imd.h"
70 #include "gromacs/listed-forces/manage-threading.h"
71 #include "gromacs/math/functions.h"
72 #include "gromacs/math/utilities.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/math/vectypes.h"
75 #include "gromacs/mdlib/checkpointhandler.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/resethandler.h"
94 #include "gromacs/mdlib/shellfc.h"
95 #include "gromacs/mdlib/sighandler.h"
96 #include "gromacs/mdlib/sim_util.h"
97 #include "gromacs/mdlib/simulationsignal.h"
98 #include "gromacs/mdlib/stophandler.h"
99 #include "gromacs/mdlib/tgroup.h"
100 #include "gromacs/mdlib/trajectory_writing.h"
101 #include "gromacs/mdlib/update.h"
102 #include "gromacs/mdlib/vcm.h"
103 #include "gromacs/mdlib/vsite.h"
104 #include "gromacs/mdtypes/awh-history.h"
105 #include "gromacs/mdtypes/awh-params.h"
106 #include "gromacs/mdtypes/commrec.h"
107 #include "gromacs/mdtypes/df_history.h"
108 #include "gromacs/mdtypes/energyhistory.h"
109 #include "gromacs/mdtypes/fcdata.h"
110 #include "gromacs/mdtypes/forcerec.h"
111 #include "gromacs/mdtypes/group.h"
112 #include "gromacs/mdtypes/inputrec.h"
113 #include "gromacs/mdtypes/interaction_const.h"
114 #include "gromacs/mdtypes/md_enums.h"
115 #include "gromacs/mdtypes/mdatom.h"
116 #include "gromacs/mdtypes/observableshistory.h"
117 #include "gromacs/mdtypes/pullhistory.h"
118 #include "gromacs/mdtypes/state.h"
119 #include "gromacs/pbcutil/mshift.h"
120 #include "gromacs/pbcutil/pbc.h"
121 #include "gromacs/pulling/output.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"
139 #include "replicaexchange.h"
142 #include "corewrap.h"
145 using gmx::SimulationSignaller;
147 void gmx::Integrator::do_md()
149 // TODO Historically, the EM and MD "integrators" used different
150 // names for the t_inputrec *parameter, but these must have the
151 // same name, now that it's a member of a struct. We use this ir
152 // alias to avoid a large ripple of nearly useless changes.
153 // t_inputrec is being replaced by IMdpOptionsProvider, so this
154 // will go away eventually.
155 t_inputrec *ir = inputrec;
156 gmx_mdoutf *outf = nullptr;
157 int64_t step, step_rel;
158 double t, t0, lam0[efptNR];
159 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
160 gmx_bool bNS, bNStList, bSimAnn, bStopCM,
161 bFirstStep, bInitStep, bLastStep = FALSE;
162 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
163 gmx_bool do_ene, do_log, do_verbose;
164 gmx_bool bMasterState;
165 int force_flags, cglo_flags;
166 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
169 matrix parrinellorahmanMu, M;
170 gmx_repl_ex_t repl_ex = nullptr;
172 t_mdebin *mdebin = nullptr;
173 gmx_enerdata_t *enerd;
174 PaddedVector<gmx::RVec> f {};
175 gmx_global_stat_t gstat;
176 gmx_update_t *upd = nullptr;
177 t_graph *graph = nullptr;
178 gmx_groups_t *groups;
179 gmx_shellfc_t *shellfc;
180 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
181 gmx_bool bTemp, bPres, bTrotter;
183 rvec *cbuf = nullptr;
190 real saved_conserved_quantity = 0;
193 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
195 /* PME load balancing data for GPU kernels */
196 pme_load_balancing_t *pme_loadbal = nullptr;
197 gmx_bool bPMETune = FALSE;
198 gmx_bool bPMETunePrinting = FALSE;
201 gmx_bool bIMDstep = FALSE;
203 /* Domain decomposition could incorrectly miss a bonded
204 interaction, but checking for that requires a global
205 communication stage, which does not otherwise happen in DD
206 code. So we do that alongside the first global energy reduction
207 after a new DD is made. These variables handle whether the
208 check happens, and the result it returns. */
209 bool shouldCheckNumberOfBondedInteractions = false;
210 int totalNumberOfBondedInteractions = -1;
212 SimulationSignals signals;
213 // Most global communnication stages don't propagate mdrun
214 // signals, and will use this object to achieve that.
215 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
217 if (!mdrunOptions.writeConfout)
219 // This is on by default, and the main known use case for
220 // turning it off is for convenience in benchmarking, which is
221 // something that should not show up in the general user
223 GMX_LOG(mdlog.info).asParagraph().
224 appendText("The -noconfout functionality is deprecated, and may be removed in a future version.");
227 /* md-vv uses averaged full step velocities for T-control
228 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
229 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
230 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
232 const bool bRerunMD = false;
233 int nstglobalcomm = mdrunOptions.globalCommunicationInterval;
235 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir, cr);
236 bGStatEveryStep = (nstglobalcomm == 1);
238 groups = &top_global->groups;
240 std::unique_ptr<EssentialDynamics> ed = nullptr;
241 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
243 /* Initialize essential dynamics sampling */
244 ed = init_edsam(mdlog,
245 opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm),
248 state_global, observablesHistory,
249 oenv, mdrunOptions.continuationOptions.appendFiles);
253 init_md(fplog, cr, outputProvider, ir, oenv, mdrunOptions,
254 &t, &t0, state_global, lam0,
255 nrnb, top_global, &upd, deform,
256 nfile, fnm, &outf, &mdebin,
257 force_vir, shake_vir, total_vir, pres, mu_tot, &bSimAnn, wcycle);
259 /* Energy terms and groups */
261 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
264 /* Kinetic energy data */
265 std::unique_ptr<gmx_ekindata_t> eKinData = compat::make_unique<gmx_ekindata_t>();
266 gmx_ekindata_t *ekind = eKinData.get();
267 init_ekindata(fplog, top_global, &(ir->opts), ekind);
268 /* Copy the cos acceleration to the groups struct */
269 ekind->cosacc.cos_accel = ir->cos_accel;
271 gstat = global_stat_init(ir);
273 /* Check for polarizable models and flexible constraints */
274 shellfc = init_shell_flexcon(fplog,
275 top_global, constr ? constr->numFlexibleConstraints() : 0,
276 ir->nstcalcenergy, DOMAINDECOMP(cr));
279 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
280 if ((io > 2000) && MASTER(cr))
283 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
288 /* Set up interactive MD (IMD) */
289 init_IMD(ir, cr, ms, top_global, fplog, ir->nstcalcenergy,
290 MASTER(cr) ? state_global->x.rvec_array() : nullptr,
291 nfile, fnm, oenv, mdrunOptions);
293 // Local state only becomes valid now.
294 std::unique_ptr<t_state> stateInstance;
297 if (DOMAINDECOMP(cr))
299 dd_init_local_top(*top_global, &top);
301 stateInstance = compat::make_unique<t_state>();
302 state = stateInstance.get();
303 dd_init_local_state(cr->dd, state_global, state);
305 /* Distribute the charge groups over the nodes from the master node */
306 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
307 state_global, *top_global, ir,
308 state, &f, mdAtoms, &top, fr,
310 nrnb, nullptr, FALSE);
311 shouldCheckNumberOfBondedInteractions = true;
312 update_realloc(upd, state->natoms);
316 state_change_natoms(state_global, state_global->natoms);
317 f.resizeWithPadding(state_global->natoms);
318 /* Copy the pointer to the global state */
319 state = state_global;
321 /* Generate and initialize new topology */
322 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr,
323 &graph, mdAtoms, constr, vsite, shellfc);
325 update_realloc(upd, state->natoms);
328 auto mdatoms = mdAtoms->mdatoms();
330 // NOTE: The global state is no longer used at this point.
331 // But state_global is still used as temporary storage space for writing
332 // the global state to file and potentially for replica exchange.
333 // (Global topology should persist.)
335 update_mdatoms(mdatoms, state->lambda[efptMASS]);
337 const ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
338 bool startingFromCheckpoint = continuationOptions.startedFromCheckpoint;
342 /* Check nstexpanded here, because the grompp check was broken */
343 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
345 gmx_fatal(FARGS, "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
347 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
352 if (startingFromCheckpoint)
354 /* Update mdebin with energy history if appending to output files */
355 if (continuationOptions.appendFiles)
357 /* If no history is available (because a checkpoint is from before
358 * it was written) make a new one later, otherwise restore it.
360 if (observablesHistory->energyHistory)
362 restore_energyhistory_from_state(mdebin, observablesHistory->energyHistory.get());
365 else if (observablesHistory->energyHistory)
367 /* We might have read an energy history from checkpoint.
368 * As we are not appending, we want to restart the statistics.
369 * Free the allocated memory and reset the counts.
371 observablesHistory->energyHistory = {};
372 /* We might have read a pull history from checkpoint.
373 * We will still want to keep the statistics, so that the files
374 * can be joined and still be meaningful.
375 * This means that observablesHistory->pullHistory
376 * should not be reset.
380 if (!observablesHistory->energyHistory)
382 observablesHistory->energyHistory = compat::make_unique<energyhistory_t>();
384 if (!observablesHistory->pullHistory)
386 observablesHistory->pullHistory = compat::make_unique<PullHistory>();
388 /* Set the initial energy history in state by updating once */
389 update_energyhistory(observablesHistory->energyHistory.get(), mdebin);
392 preparePrevStepPullCom(ir, mdatoms, state, state_global, cr, startingFromCheckpoint);
394 // TODO: Remove this by converting AWH into a ForceProvider
395 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms, startingFromCheckpoint,
397 opt2fn("-awh", nfile, fnm), ir->pull_work);
399 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
400 if (useReplicaExchange && MASTER(cr))
402 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir,
405 /* PME tuning is only supported in the Verlet scheme, with PME for
406 * Coulomb. It is not supported with only LJ PME. */
407 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) &&
408 !mdrunOptions.reproducible && ir->cutoff_scheme != ecutsGROUP);
411 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box,
412 *fr->ic, *fr->nbv->listParams, fr->pmedata, use_GPU(fr->nbv),
416 if (!ir->bContinuation)
418 if (state->flags & (1 << estV))
420 auto v = makeArrayRef(state->v);
421 /* Set the velocities of vsites, shells and frozen atoms to zero */
422 for (i = 0; i < mdatoms->homenr; i++)
424 if (mdatoms->ptype[i] == eptVSite ||
425 mdatoms->ptype[i] == eptShell)
429 else if (mdatoms->cFREEZE)
431 for (m = 0; m < DIM; m++)
433 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
444 /* Constrain the initial coordinates and velocities */
445 do_constrain_first(fplog, constr, ir, mdatoms, state);
449 /* Construct the virtual sites for the initial configuration */
450 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, nullptr,
451 top.idef.iparams, top.idef.il,
452 fr->ePBC, fr->bMolPBC, cr, state->box);
456 if (ir->efep != efepNO)
458 /* Set free energy calculation frequency as the greatest common
459 * denominator of nstdhdl and repl_ex_nst. */
460 nstfep = ir->fepvals->nstdhdl;
463 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
465 if (useReplicaExchange)
467 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
471 /* Be REALLY careful about what flags you set here. You CANNOT assume
472 * this is the first step, since we might be restarting from a checkpoint,
473 * and in that case we should not do any modifications to the state.
475 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
477 if (continuationOptions.haveReadEkin)
479 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
482 cglo_flags = (CGLO_INITIALIZATION | CGLO_TEMPERATURE | CGLO_GSTAT
483 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
484 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
485 | (continuationOptions.haveReadEkin ? CGLO_READEKIN : 0));
487 bSumEkinhOld = FALSE;
489 t_vcm vcm(top_global->groups, *ir);
490 reportComRemovalInfo(fplog, vcm);
492 /* To minimize communication, compute_globals computes the COM velocity
493 * and the kinetic energy for the velocities without COM motion removed.
494 * Thus to get the kinetic energy without the COM contribution, we need
495 * to call compute_globals twice.
497 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
499 int cglo_flags_iteration = cglo_flags;
500 if (bStopCM && cgloIteration == 0)
502 cglo_flags_iteration |= CGLO_STOPCM;
503 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
505 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
506 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
507 constr, &nullSignaller, state->box,
508 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags_iteration
509 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
511 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
512 top_global, &top, state,
513 &shouldCheckNumberOfBondedInteractions);
514 if (ir->eI == eiVVAK)
516 /* a second call to get the half step temperature initialized as well */
517 /* we do the same call as above, but turn the pressure off -- internally to
518 compute_globals, this is recognized as a velocity verlet half-step
519 kinetic energy calculation. This minimized excess variables, but
520 perhaps loses some logic?*/
522 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
523 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
524 constr, &nullSignaller, state->box,
525 nullptr, &bSumEkinhOld,
526 cglo_flags & ~CGLO_PRESSURE);
529 /* Calculate the initial half step temperature, and save the ekinh_old */
530 if (!continuationOptions.startedFromCheckpoint)
532 for (i = 0; (i < ir->opts.ngtc); i++)
534 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
538 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
539 temperature control */
540 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
544 if (!ir->bContinuation)
546 if (constr && ir->eConstrAlg == econtLINCS)
549 "RMS relative constraint deviation after constraining: %.2e\n",
552 if (EI_STATE_VELOCITY(ir->eI))
554 real temp = enerd->term[F_TEMP];
557 /* Result of Ekin averaged over velocities of -half
558 * and +half step, while we only have -half step here.
562 fprintf(fplog, "Initial temperature: %g K\n", temp);
567 fprintf(stderr, "starting mdrun '%s'\n",
568 *(top_global->name));
571 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
575 sprintf(tbuf, "%s", "infinite");
577 if (ir->init_step > 0)
579 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
580 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
581 gmx_step_str(ir->init_step, sbuf2),
582 ir->init_step*ir->delta_t);
586 fprintf(stderr, "%s steps, %s ps.\n",
587 gmx_step_str(ir->nsteps, sbuf), tbuf);
589 fprintf(fplog, "\n");
592 walltime_accounting_start_time(walltime_accounting);
593 wallcycle_start(wcycle, ewcRUN);
594 print_start(fplog, cr, walltime_accounting, "mdrun");
597 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
598 int chkpt_ret = fcCheckPointParallel( cr->nodeid,
602 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
606 /***********************************************************
610 ************************************************************/
613 /* Skip the first Nose-Hoover integration when we get the state from tpx */
614 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
615 bSumEkinhOld = FALSE;
617 bNeedRepartition = FALSE;
619 bool simulationsShareState = false;
620 int nstSignalComm = nstglobalcomm;
622 // TODO This implementation of ensemble orientation restraints is nasty because
623 // a user can't just do multi-sim with single-sim orientation restraints.
624 bool usingEnsembleRestraints = (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
625 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
627 // Replica exchange, ensemble restraints and AWH need all
628 // simulations to remain synchronized, so they need
629 // checkpoints and stop conditions to act on the same step, so
630 // the propagation of such signals must take place between
631 // simulations, not just within simulations.
632 // TODO: Make algorithm initializers set these flags.
633 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
635 if (simulationsShareState)
637 // Inter-simulation signal communication does not need to happen
638 // often, so we use a minimum of 200 steps to reduce overhead.
639 const int c_minimumInterSimulationSignallingInterval = 200;
640 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1)/nstglobalcomm)*nstglobalcomm;
644 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
645 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
646 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
647 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
649 auto checkpointHandler = compat::make_unique<CheckpointHandler>(
650 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]),
651 simulationsShareState, ir->nstlist == 0, MASTER(cr),
652 mdrunOptions.writeConfout, mdrunOptions.checkpointOptions.period);
654 const bool resetCountersIsLocal = true;
655 auto resetHandler = compat::make_unique<ResetHandler>(
656 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]), !resetCountersIsLocal,
657 ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
658 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
660 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion = (DOMAINDECOMP(cr) ? DdOpenBalanceRegionBeforeForceComputation::yes : DdOpenBalanceRegionBeforeForceComputation::no);
661 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion = (DOMAINDECOMP(cr) ? DdCloseBalanceRegionAfterForceComputation::yes : DdCloseBalanceRegionAfterForceComputation::no);
663 step = ir->init_step;
666 // TODO extract this to new multi-simulation module
667 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
669 if (!multisim_int_all_are_equal(ms, ir->nsteps))
671 GMX_LOG(mdlog.warning).appendText(
672 "Note: The number of steps is not consistent across multi simulations,\n"
673 "but we are proceeding anyway!");
675 if (!multisim_int_all_are_equal(ms, ir->init_step))
677 GMX_LOG(mdlog.warning).appendText(
678 "Note: The initial step is not consistent across multi simulations,\n"
679 "but we are proceeding anyway!");
683 /* and stop now if we should */
684 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
688 /* Determine if this is a neighbor search step */
689 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
691 if (bPMETune && bNStList)
693 /* PME grid + cut-off optimization with GPUs or PME nodes */
694 pme_loadbal_do(pme_loadbal, cr,
695 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
703 wallcycle_start(wcycle, ewcSTEP);
705 bLastStep = (step_rel == ir->nsteps);
706 t = t0 + step*ir->delta_t;
708 // TODO Refactor this, so that nstfep does not need a default value of zero
709 if (ir->efep != efepNO || ir->bSimTemp)
711 /* find and set the current lambdas */
712 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
714 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
715 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
716 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
717 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
720 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep &&
721 do_per_step(step, replExParams.exchangeInterval));
725 update_annealing_target_temp(ir, t, upd);
728 /* Stop Center of Mass motion */
729 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
731 /* Determine whether or not to do Neighbour Searching */
732 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
734 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
736 /* do_log triggers energy and virial calculation. Because this leads
737 * to different code paths, forces can be different. Thus for exact
738 * continuation we should avoid extra log output.
739 * Note that the || bLastStep can result in non-exact continuation
740 * beyond the last step. But we don't consider that to be an issue.
742 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep;
743 do_verbose = mdrunOptions.verbose &&
744 (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
746 if (bNS && !(bFirstStep && ir->bContinuation))
748 bMasterState = FALSE;
749 /* Correct the new box if it is too skewed */
750 if (inputrecDynamicBox(ir))
752 if (correct_box(fplog, step, state->box, graph))
757 if (DOMAINDECOMP(cr) && bMasterState)
759 dd_collect_state(cr->dd, state, state_global);
762 if (DOMAINDECOMP(cr))
764 /* Repartition the domain decomposition */
765 dd_partition_system(fplog, mdlog, step, cr,
766 bMasterState, nstglobalcomm,
767 state_global, *top_global, ir,
768 state, &f, mdAtoms, &top, fr,
771 do_verbose && !bPMETunePrinting);
772 shouldCheckNumberOfBondedInteractions = true;
773 update_realloc(upd, state->natoms);
777 if (MASTER(cr) && do_log)
779 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
782 if (ir->efep != efepNO)
784 update_mdatoms(mdatoms, state->lambda[efptMASS]);
790 /* We need the kinetic energy at minus the half step for determining
791 * the full step kinetic energy and possibly for T-coupling.*/
792 /* This may not be quite working correctly yet . . . . */
793 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
794 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
795 constr, &nullSignaller, state->box,
796 &totalNumberOfBondedInteractions, &bSumEkinhOld,
797 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
798 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
799 top_global, &top, state,
800 &shouldCheckNumberOfBondedInteractions);
802 clear_mat(force_vir);
804 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
806 /* Determine the energy and pressure:
807 * at nstcalcenergy steps and at energy output steps (set below).
809 if (EI_VV(ir->eI) && (!bInitStep))
811 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
812 bCalcVir = bCalcEnerStep ||
813 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
817 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
818 bCalcVir = bCalcEnerStep ||
819 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
821 bCalcEner = bCalcEnerStep;
823 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
825 if (do_ene || do_log || bDoReplEx)
831 /* Do we need global communication ? */
832 bGStat = (bCalcVir || bCalcEner || bStopCM ||
833 do_per_step(step, nstglobalcomm) ||
834 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
836 force_flags = (GMX_FORCE_STATECHANGED |
837 ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0) |
838 GMX_FORCE_ALLFORCES |
839 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
840 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
841 (bDoFEP ? GMX_FORCE_DHDL : 0)
846 /* Now is the time to relax the shells */
847 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
848 enforcedRotation, step,
849 ir, bNS, force_flags, &top,
851 state, f.arrayRefWithPadding(), force_vir, mdatoms,
852 nrnb, wcycle, graph, groups,
853 shellfc, fr, ppForceWorkload, t, mu_tot,
855 ddOpenBalanceRegion, ddCloseBalanceRegion);
859 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
860 (or the AWH update will be performed twice for one step when continuing). It would be best to
861 call this update function from do_md_trajectory_writing but that would occur after do_force.
862 One would have to divide the update_awh function into one function applying the AWH force
863 and one doing the AWH bias update. The update AWH bias function could then be called after
864 do_md_trajectory_writing (then containing update_awh_history).
865 The checkpointing will in the future probably moved to the start of the md loop which will
866 rid of this issue. */
867 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
869 awh->updateHistory(state_global->awhHistory.get());
872 /* The coordinates (x) are shifted (to get whole molecules)
874 * This is parallellized as well, and does communication too.
875 * Check comments in sim_util.c
877 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation,
878 step, nrnb, wcycle, &top, groups,
879 state->box, state->x.arrayRefWithPadding(), &state->hist,
880 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
881 state->lambda, graph,
882 fr, ppForceWorkload, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
883 (bNS ? GMX_FORCE_NS : 0) | force_flags,
884 ddOpenBalanceRegion, ddCloseBalanceRegion);
887 if (EI_VV(ir->eI) && !startingFromCheckpoint)
888 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
890 rvec *vbuf = nullptr;
892 wallcycle_start(wcycle, ewcUPDATE);
893 if (ir->eI == eiVV && bInitStep)
895 /* if using velocity verlet with full time step Ekin,
896 * take the first half step only to compute the
897 * virial for the first step. From there,
898 * revert back to the initial coordinates
899 * so that the input is actually the initial step.
901 snew(vbuf, state->natoms);
902 copy_rvecn(state->v.rvec_array(), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
906 /* this is for NHC in the Ekin(t+dt/2) version of vv */
907 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
910 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
911 ekind, M, upd, etrtVELOCITY1,
914 wallcycle_stop(wcycle, ewcUPDATE);
915 constrain_velocities(step, nullptr,
919 bCalcVir, do_log, do_ene);
920 wallcycle_start(wcycle, ewcUPDATE);
921 /* if VV, compute the pressure and constraints */
922 /* For VV2, we strictly only need this if using pressure
923 * control, but we really would like to have accurate pressures
925 * Think about ways around this in the future?
926 * For now, keep this choice in comments.
928 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
929 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
931 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
932 if (bCalcEner && ir->eI == eiVVAK)
936 /* for vv, the first half of the integration actually corresponds to the previous step.
937 So we need information from the last step in the first half of the integration */
938 if (bGStat || do_per_step(step-1, nstglobalcomm))
940 wallcycle_stop(wcycle, ewcUPDATE);
941 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
942 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
943 constr, &nullSignaller, state->box,
944 &totalNumberOfBondedInteractions, &bSumEkinhOld,
945 (bGStat ? CGLO_GSTAT : 0)
946 | (bCalcEner ? CGLO_ENERGY : 0)
947 | (bTemp ? CGLO_TEMPERATURE : 0)
948 | (bPres ? CGLO_PRESSURE : 0)
949 | (bPres ? CGLO_CONSTRAINT : 0)
950 | (bStopCM ? CGLO_STOPCM : 0)
951 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
954 /* explanation of above:
955 a) We compute Ekin at the full time step
956 if 1) we are using the AveVel Ekin, and it's not the
957 initial step, or 2) if we are using AveEkin, but need the full
958 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
959 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
960 EkinAveVel because it's needed for the pressure */
961 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
962 top_global, &top, state,
963 &shouldCheckNumberOfBondedInteractions);
964 wallcycle_start(wcycle, ewcUPDATE);
966 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
971 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
972 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
974 /* TODO This is only needed when we're about to write
975 * a checkpoint, because we use it after the restart
976 * (in a kludge?). But what should we be doing if
977 * startingFromCheckpoint or bInitStep are true? */
978 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
980 copy_mat(shake_vir, state->svir_prev);
981 copy_mat(force_vir, state->fvir_prev);
983 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
985 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
986 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
987 enerd->term[F_EKIN] = trace(ekind->ekin);
992 wallcycle_stop(wcycle, ewcUPDATE);
993 /* We need the kinetic energy at minus the half step for determining
994 * the full step kinetic energy and possibly for T-coupling.*/
995 /* This may not be quite working correctly yet . . . . */
996 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
997 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
998 constr, &nullSignaller, state->box,
999 nullptr, &bSumEkinhOld,
1000 CGLO_GSTAT | CGLO_TEMPERATURE);
1001 wallcycle_start(wcycle, ewcUPDATE);
1004 /* if it's the initial step, we performed this first step just to get the constraint virial */
1005 if (ir->eI == eiVV && bInitStep)
1007 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1010 wallcycle_stop(wcycle, ewcUPDATE);
1013 /* compute the conserved quantity */
1016 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1019 last_ekin = enerd->term[F_EKIN];
1021 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1023 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1025 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1026 if (ir->efep != efepNO)
1028 sum_dhdl(enerd, state->lambda, ir->fepvals);
1032 /* ######## END FIRST UPDATE STEP ############## */
1033 /* ######## If doing VV, we now have v(dt) ###### */
1036 /* perform extended ensemble sampling in lambda - we don't
1037 actually move to the new state before outputting
1038 statistics, but if performing simulated tempering, we
1039 do update the velocities and the tau_t. */
1041 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, state->v.rvec_array(), mdatoms);
1042 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1045 copy_df_history(state_global->dfhist, state->dfhist);
1049 /* Now we have the energies and forces corresponding to the
1050 * coordinates at time t. We must output all of this before
1053 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1054 ir, state, state_global, observablesHistory,
1056 outf, mdebin, ekind, f,
1057 checkpointHandler->isCheckpointingStep(),
1058 bRerunMD, bLastStep,
1059 mdrunOptions.writeConfout,
1061 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1062 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x.rvec_array(), ir, t, wcycle);
1064 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1065 if (startingFromCheckpoint && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1067 copy_mat(state->svir_prev, shake_vir);
1068 copy_mat(state->fvir_prev, force_vir);
1071 stopHandler->setSignal();
1072 resetHandler->setSignal(walltime_accounting);
1074 if (bGStat || !PAR(cr))
1076 /* In parallel we only have to check for checkpointing in steps
1077 * where we do global communication,
1078 * otherwise the other nodes don't know.
1080 checkpointHandler->setSignal(walltime_accounting);
1083 /* ######### START SECOND UPDATE STEP ################# */
1085 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1088 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1090 gmx_bool bIfRandomize;
1091 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, upd, constr);
1092 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1093 if (constr && bIfRandomize)
1095 constrain_velocities(step, nullptr,
1099 bCalcVir, do_log, do_ene);
1102 /* Box is changed in update() when we do pressure coupling,
1103 * but we should still use the old box for energy corrections and when
1104 * writing it to the energy file, so it matches the trajectory files for
1105 * the same timestep above. Make a copy in a separate array.
1107 copy_mat(state->box, lastbox);
1111 wallcycle_start(wcycle, ewcUPDATE);
1112 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1115 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1116 /* We can only do Berendsen coupling after we have summed
1117 * the kinetic energy or virial. Since the happens
1118 * in global_state after update, we should only do it at
1119 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1124 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1125 update_pcouple_before_coordinates(fplog, step, ir, state,
1126 parrinellorahmanMu, M,
1132 /* velocity half-step update */
1133 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
1134 ekind, M, upd, etrtVELOCITY2,
1138 /* Above, initialize just copies ekinh into ekin,
1139 * it doesn't copy position (for VV),
1140 * and entire integrator for MD.
1143 if (ir->eI == eiVVAK)
1145 /* We probably only need md->homenr, not state->natoms */
1146 if (state->natoms > cbuf_nalloc)
1148 cbuf_nalloc = state->natoms;
1149 srenew(cbuf, cbuf_nalloc);
1151 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1154 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
1155 ekind, M, upd, etrtPOSITION, cr, constr);
1156 wallcycle_stop(wcycle, ewcUPDATE);
1158 constrain_coordinates(step, &dvdl_constr, state,
1161 bCalcVir, do_log, do_ene);
1162 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state,
1163 cr, nrnb, wcycle, upd, constr, do_log, do_ene);
1164 finish_update(ir, mdatoms,
1166 nrnb, wcycle, upd, constr);
1168 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1170 updatePrevStepPullCom(ir->pull_work, state);
1173 if (ir->eI == eiVVAK)
1175 /* erase F_EKIN and F_TEMP here? */
1176 /* just compute the kinetic energy at the half step to perform a trotter step */
1177 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
1178 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1179 constr, &nullSignaller, lastbox,
1180 nullptr, &bSumEkinhOld,
1181 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1183 wallcycle_start(wcycle, ewcUPDATE);
1184 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1185 /* now we know the scaling, we can compute the positions again again */
1186 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1188 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
1189 ekind, M, upd, etrtPOSITION, cr, constr);
1190 wallcycle_stop(wcycle, ewcUPDATE);
1192 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1193 /* are the small terms in the shake_vir here due
1194 * to numerical errors, or are they important
1195 * physically? I'm thinking they are just errors, but not completely sure.
1196 * For now, will call without actually constraining, constr=NULL*/
1197 finish_update(ir, mdatoms,
1199 nrnb, wcycle, upd, nullptr);
1203 /* this factor or 2 correction is necessary
1204 because half of the constraint force is removed
1205 in the vv step, so we have to double it. See
1206 the Redmine issue #1255. It is not yet clear
1207 if the factor of 2 is exact, or just a very
1208 good approximation, and this will be
1209 investigated. The next step is to see if this
1210 can be done adding a dhdl contribution from the
1211 rattle step, but this is somewhat more
1212 complicated with the current code. Will be
1213 investigated, hopefully for 4.6.3. However,
1214 this current solution is much better than
1215 having it completely wrong.
1217 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1221 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1224 if (vsite != nullptr)
1226 wallcycle_start(wcycle, ewcVSITECONSTR);
1227 if (graph != nullptr)
1229 shift_self(graph, state->box, state->x.rvec_array());
1231 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, state->v.rvec_array(),
1232 top.idef.iparams, top.idef.il,
1233 fr->ePBC, fr->bMolPBC, cr, state->box);
1235 if (graph != nullptr)
1237 unshift_self(graph, state->box, state->x.rvec_array());
1239 wallcycle_stop(wcycle, ewcVSITECONSTR);
1242 /* ############## IF NOT VV, Calculate globals HERE ############ */
1243 /* With Leap-Frog we can skip compute_globals at
1244 * non-communication steps, but we need to calculate
1245 * the kinetic energy one step before communication.
1248 // Organize to do inter-simulation signalling on steps if
1249 // and when algorithms require it.
1250 bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1252 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1254 // Since we're already communicating at this step, we
1255 // can propagate intra-simulation signals. Note that
1256 // check_nstglobalcomm has the responsibility for
1257 // choosing the value of nstglobalcomm that is one way
1258 // bGStat becomes true, so we can't get into a
1259 // situation where e.g. checkpointing can't be
1261 bool doIntraSimSignal = true;
1262 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1264 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
1265 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1268 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1269 (bGStat ? CGLO_GSTAT : 0)
1270 | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1271 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1272 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1273 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
1275 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1277 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1278 top_global, &top, state,
1279 &shouldCheckNumberOfBondedInteractions);
1283 /* ############# END CALC EKIN AND PRESSURE ################# */
1285 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1286 the virial that should probably be addressed eventually. state->veta has better properies,
1287 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1288 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1290 if (ir->efep != efepNO && !EI_VV(ir->eI))
1292 /* Sum up the foreign energy and dhdl terms for md and sd.
1293 Currently done every step so that dhdl is correct in the .edr */
1294 sum_dhdl(enerd, state->lambda, ir->fepvals);
1297 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1298 pres, force_vir, shake_vir,
1302 /* ################# END UPDATE STEP 2 ################# */
1303 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1305 /* The coordinates (x) were unshifted in update */
1308 /* We will not sum ekinh_old,
1309 * so signal that we still have to do it.
1311 bSumEkinhOld = TRUE;
1316 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1318 /* use the directly determined last velocity, not actually the averaged half steps */
1319 if (bTrotter && ir->eI == eiVV)
1321 enerd->term[F_EKIN] = last_ekin;
1323 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1325 if (integratorHasConservedEnergyQuantity(ir))
1329 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1333 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1336 /* ######### END PREPARING EDR OUTPUT ########### */
1342 if (fplog && do_log && bDoExpanded)
1344 /* only needed if doing expanded ensemble */
1345 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1346 state_global->dfhist, state->fep_state, ir->nstlog, step);
1350 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1351 t, mdatoms->tmass, enerd, state,
1352 ir->fepvals, ir->expandedvals, lastbox,
1353 shake_vir, force_vir, total_vir, pres,
1354 ekind, mu_tot, constr);
1358 upd_mdebin_step(mdebin);
1361 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1362 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1364 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1366 eprNORMAL, mdebin, fcd, groups, &(ir->opts), awh.get());
1370 pull_print_output(ir->pull_work, step, t);
1373 if (do_per_step(step, ir->nstlog))
1375 if (fflush(fplog) != 0)
1377 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1383 /* Have to do this part _after_ outputting the logfile and the edr file */
1384 /* Gets written into the state at the beginning of next loop*/
1385 state->fep_state = lamnew;
1387 /* Print the remaining wall clock time for the run */
1388 if (isMasterSimMasterRank(ms, cr) &&
1389 (do_verbose || gmx_got_usr_signal()) &&
1394 fprintf(stderr, "\n");
1396 print_time(stderr, walltime_accounting, step, ir, cr);
1399 /* Ion/water position swapping.
1400 * Not done in last step since trajectory writing happens before this call
1401 * in the MD loop and exchanges would be lost anyway. */
1402 bNeedRepartition = FALSE;
1403 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1404 do_per_step(step, ir->swap->nstswap))
1406 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1407 as_rvec_array(state->x.data()),
1409 MASTER(cr) && mdrunOptions.verbose,
1412 if (bNeedRepartition && DOMAINDECOMP(cr))
1414 dd_collect_state(cr->dd, state, state_global);
1418 /* Replica exchange */
1422 bExchanged = replica_exchange(fplog, cr, ms, repl_ex,
1423 state_global, enerd,
1427 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1429 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1,
1430 state_global, *top_global, ir,
1431 state, &f, mdAtoms, &top, fr,
1433 nrnb, wcycle, FALSE);
1434 shouldCheckNumberOfBondedInteractions = true;
1435 update_realloc(upd, state->natoms);
1440 startingFromCheckpoint = false;
1442 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1443 /* With all integrators, except VV, we need to retain the pressure
1444 * at the current step for coupling at the next step.
1446 if ((state->flags & (1<<estPRES_PREV)) &&
1448 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1450 /* Store the pressure in t_state for pressure coupling
1451 * at the next MD step.
1453 copy_mat(pres, state->pres_prev);
1456 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1458 if ( (membed != nullptr) && (!bLastStep) )
1460 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1463 cycles = wallcycle_stop(wcycle, ewcSTEP);
1464 if (DOMAINDECOMP(cr) && wcycle)
1466 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1469 /* increase the MD step number */
1473 resetHandler->resetCounters(
1474 step, step_rel, mdlog, fplog, cr, (use_GPU(fr->nbv) ? fr->nbv : nullptr),
1475 nrnb, fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1477 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1478 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1481 /* End of main MD loop */
1483 /* Closing TNG files can include compressing data. Therefore it is good to do that
1484 * before stopping the time measurements. */
1485 mdoutf_tng_close(outf);
1487 /* Stop measuring walltime */
1488 walltime_accounting_end_time(walltime_accounting);
1490 if (!thisRankHasDuty(cr, DUTY_PME))
1492 /* Tell the PME only node to finish */
1493 gmx_pme_send_finish(cr);
1498 if (ir->nstcalcenergy > 0)
1500 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1501 eprAVER, mdebin, fcd, groups, &(ir->opts), awh.get());
1504 done_mdebin(mdebin);
1509 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
1512 done_shellfc(fplog, shellfc, step_rel);
1514 if (useReplicaExchange && MASTER(cr))
1516 print_replica_exchange_statistics(fplog, repl_ex);
1519 // Clean up swapcoords
1520 if (ir->eSwapCoords != eswapNO)
1522 finish_swapcoords(ir->swap);
1525 /* IMD cleanup, if bIMD is TRUE. */
1526 IMD_finalize(ir->bIMD, ir->imd);
1528 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1530 destroy_enerdata(enerd);
1533 global_stat_destroy(gstat);