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/domdec/collect.h"
57 #include "gromacs/domdec/dlbtiming.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/mdsetup.h"
62 #include "gromacs/domdec/partition.h"
63 #include "gromacs/essentialdynamics/edsam.h"
64 #include "gromacs/ewald/pme.h"
65 #include "gromacs/ewald/pme_load_balancing.h"
66 #include "gromacs/fileio/trxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/gpu_utils/gpu_utils.h"
70 #include "gromacs/imd/imd.h"
71 #include "gromacs/listed_forces/manage_threading.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/utilities.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/math/vectypes.h"
76 #include "gromacs/mdlib/checkpointhandler.h"
77 #include "gromacs/mdlib/compute_io.h"
78 #include "gromacs/mdlib/constr.h"
79 #include "gromacs/mdlib/ebin.h"
80 #include "gromacs/mdlib/enerdata_utils.h"
81 #include "gromacs/mdlib/energyoutput.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/mdoutf.h"
89 #include "gromacs/mdlib/membed.h"
90 #include "gromacs/mdlib/resethandler.h"
91 #include "gromacs/mdlib/sighandler.h"
92 #include "gromacs/mdlib/simulationsignal.h"
93 #include "gromacs/mdlib/stat.h"
94 #include "gromacs/mdlib/stophandler.h"
95 #include "gromacs/mdlib/tgroup.h"
96 #include "gromacs/mdlib/trajectory_writing.h"
97 #include "gromacs/mdlib/update.h"
98 #include "gromacs/mdlib/update_constrain_cuda.h"
99 #include "gromacs/mdlib/vcm.h"
100 #include "gromacs/mdlib/vsite.h"
101 #include "gromacs/mdrunutility/handlerestart.h"
102 #include "gromacs/mdrunutility/multisim.h"
103 #include "gromacs/mdrunutility/printtime.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/mdrunoptions.h"
117 #include "gromacs/mdtypes/observableshistory.h"
118 #include "gromacs/mdtypes/pullhistory.h"
119 #include "gromacs/mdtypes/state.h"
120 #include "gromacs/nbnxm/nbnxm.h"
121 #include "gromacs/pbcutil/mshift.h"
122 #include "gromacs/pbcutil/pbc.h"
123 #include "gromacs/pulling/output.h"
124 #include "gromacs/pulling/pull.h"
125 #include "gromacs/swap/swapcoords.h"
126 #include "gromacs/timing/wallcycle.h"
127 #include "gromacs/timing/walltime_accounting.h"
128 #include "gromacs/topology/atoms.h"
129 #include "gromacs/topology/idef.h"
130 #include "gromacs/topology/mtop_util.h"
131 #include "gromacs/topology/topology.h"
132 #include "gromacs/trajectory/trajectoryframe.h"
133 #include "gromacs/utility/basedefinitions.h"
134 #include "gromacs/utility/cstringutil.h"
135 #include "gromacs/utility/fatalerror.h"
136 #include "gromacs/utility/logger.h"
137 #include "gromacs/utility/real.h"
138 #include "gromacs/utility/smalloc.h"
140 #include "replicaexchange.h"
142 #include "simulator.h"
145 #include "corewrap.h"
148 using gmx::SimulationSignaller;
150 //! Whether the GPU versions of Leap-Frog integrator and LINCS and SHAKE constraints
151 static const bool c_useGpuUpdateConstrain = (getenv("GMX_UPDATE_CONSTRAIN_GPU") != nullptr);
153 void gmx::Simulator::do_md()
155 // TODO Historically, the EM and MD "integrators" used different
156 // names for the t_inputrec *parameter, but these must have the
157 // same name, now that it's a member of a struct. We use this ir
158 // alias to avoid a large ripple of nearly useless changes.
159 // t_inputrec is being replaced by IMdpOptionsProvider, so this
160 // will go away eventually.
161 t_inputrec *ir = inputrec;
162 int64_t step, step_rel;
163 double t, t0 = ir->init_t, lam0[efptNR];
164 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
165 gmx_bool bNS, bNStList, bStopCM,
166 bFirstStep, bInitStep, bLastStep = FALSE;
167 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
168 gmx_bool do_ene, do_log, do_verbose;
169 gmx_bool bMasterState;
170 int force_flags, cglo_flags;
171 tensor force_vir = {{0}}, shake_vir = {{0}}, total_vir = {{0}},
172 tmp_vir = {{0}}, pres = {{0}};
175 matrix parrinellorahmanMu, M;
176 gmx_repl_ex_t repl_ex = nullptr;
178 PaddedVector<gmx::RVec> f {};
179 gmx_global_stat_t gstat;
180 t_graph *graph = nullptr;
181 gmx_shellfc_t *shellfc;
182 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
183 gmx_bool bTemp, bPres, bTrotter;
185 rvec *cbuf = nullptr;
192 real saved_conserved_quantity = 0;
195 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
197 /* PME load balancing data for GPU kernels */
198 gmx_bool bPMETune = FALSE;
199 gmx_bool bPMETunePrinting = FALSE;
201 bool bInteractiveMDstep = 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;
234 int nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
235 bGStatEveryStep = (nstglobalcomm == 1);
237 SimulationGroups *groups = &top_global->groups;
239 std::unique_ptr<EssentialDynamics> ed = nullptr;
240 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
242 /* Initialize essential dynamics sampling */
243 ed = init_edsam(mdlog,
244 opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm),
247 state_global, observablesHistory,
252 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
253 Update upd(ir, deform);
254 bool doSimulatedAnnealing = initSimulatedAnnealing(ir, &upd);
255 if (startingBehavior != StartingBehavior::RestartWithAppending)
257 pleaseCiteCouplingAlgorithms(fplog, *ir);
260 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, oenv, wcycle,
261 StartingBehavior::NewSimulation);
262 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work, mdoutf_get_fp_dhdl(outf), false);
264 /* Kinetic energy data */
265 std::unique_ptr<gmx_ekindata_t> eKinData = std::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, energyOutput.numEnergyTerms(), 1);
280 if ((io > 2000) && MASTER(cr))
283 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
288 // Local state only becomes valid now.
289 std::unique_ptr<t_state> stateInstance;
293 auto mdatoms = mdAtoms->mdatoms();
295 std::unique_ptr<UpdateConstrainCuda> integrator;
297 if (DOMAINDECOMP(cr))
299 dd_init_local_top(*top_global, &top);
301 stateInstance = std::make_unique<t_state>();
302 state = stateInstance.get();
303 if (fr->nbv->useGpu())
305 changePinningPolicy(&state->x, gmx::PinningPolicy::PinnedIfSupported);
307 dd_init_local_state(cr->dd, state_global, state);
309 /* Distribute the charge groups over the nodes from the master node */
310 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
311 state_global, *top_global, ir, imdSession,
313 state, &f, mdAtoms, &top, fr,
315 nrnb, nullptr, FALSE);
316 shouldCheckNumberOfBondedInteractions = true;
317 upd.setNumAtoms(state->natoms);
321 state_change_natoms(state_global, state_global->natoms);
322 f.resizeWithPadding(state_global->natoms);
323 /* Copy the pointer to the global state */
324 state = state_global;
326 /* Generate and initialize new topology */
327 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr,
328 &graph, mdAtoms, constr, vsite, shellfc);
330 upd.setNumAtoms(state->natoms);
332 if (c_useGpuUpdateConstrain)
334 GMX_RELEASE_ASSERT(ir->eI == eiMD, "Only md integrator is supported on the GPU.");
335 GMX_RELEASE_ASSERT(ir->etc == etcNO, "Temperature coupling is not supported on the GPU.");
336 GMX_RELEASE_ASSERT(ir->epc == epcNO, "Pressure coupling is not supported on the GPU.");
337 GMX_RELEASE_ASSERT(!mdatoms->haveVsites, "Virtual sites are not supported on the GPU");
338 GMX_LOG(mdlog.info).asParagraph().
339 appendText("Using CUDA GPU-based update and constraints module.");
340 integrator = std::make_unique<UpdateConstrainCuda>(state->natoms, *ir, *top_global);
341 integrator->set(top.idef, *mdatoms);
343 set_pbc(&pbc, epbcXYZ, state->box);
344 integrator->setPbc(&pbc);
349 // NOTE: The global state is no longer used at this point.
350 // But state_global is still used as temporary storage space for writing
351 // the global state to file and potentially for replica exchange.
352 // (Global topology should persist.)
354 update_mdatoms(mdatoms, state->lambda[efptMASS]);
358 /* Check nstexpanded here, because the grompp check was broken */
359 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
361 gmx_fatal(FARGS, "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
363 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation,
369 if (startingBehavior != StartingBehavior::NewSimulation)
371 /* Restore from energy history if appending to output files */
372 if (startingBehavior == StartingBehavior::RestartWithAppending)
374 /* If no history is available (because a checkpoint is from before
375 * it was written) make a new one later, otherwise restore it.
377 if (observablesHistory->energyHistory)
379 energyOutput.restoreFromEnergyHistory(*observablesHistory->energyHistory);
382 else if (observablesHistory->energyHistory)
384 /* We might have read an energy history from checkpoint.
385 * As we are not appending, we want to restart the statistics.
386 * Free the allocated memory and reset the counts.
388 observablesHistory->energyHistory = {};
389 /* We might have read a pull history from checkpoint.
390 * We will still want to keep the statistics, so that the files
391 * can be joined and still be meaningful.
392 * This means that observablesHistory->pullHistory
393 * should not be reset.
397 if (!observablesHistory->energyHistory)
399 observablesHistory->energyHistory = std::make_unique<energyhistory_t>();
401 if (!observablesHistory->pullHistory)
403 observablesHistory->pullHistory = std::make_unique<PullHistory>();
405 /* Set the initial energy history */
406 energyOutput.fillEnergyHistory(observablesHistory->energyHistory.get());
409 preparePrevStepPullCom(ir, pull_work, mdatoms, state, state_global, cr,
410 startingBehavior != StartingBehavior::NewSimulation);
412 // TODO: Remove this by converting AWH into a ForceProvider
413 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms, startingBehavior != StartingBehavior::NewSimulation,
415 opt2fn("-awh", nfile, fnm), pull_work);
417 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
418 if (useReplicaExchange && MASTER(cr))
420 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir,
423 /* PME tuning is only supported in the Verlet scheme, with PME for
424 * Coulomb. It is not supported with only LJ PME. */
425 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) &&
426 !mdrunOptions.reproducible && ir->cutoff_scheme != ecutsGROUP);
428 pme_load_balancing_t *pme_loadbal = nullptr;
431 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box,
432 *fr->ic, *fr->nbv, fr->pmedata, fr->nbv->useGpu(),
436 if (!ir->bContinuation)
438 if (state->flags & (1 << estV))
440 auto v = makeArrayRef(state->v);
441 /* Set the velocities of vsites, shells and frozen atoms to zero */
442 for (i = 0; i < mdatoms->homenr; i++)
444 if (mdatoms->ptype[i] == eptVSite ||
445 mdatoms->ptype[i] == eptShell)
449 else if (mdatoms->cFREEZE)
451 for (m = 0; m < DIM; m++)
453 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
464 /* Constrain the initial coordinates and velocities */
465 do_constrain_first(fplog, constr, ir, mdatoms, state);
469 /* Construct the virtual sites for the initial configuration */
470 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, nullptr,
471 top.idef.iparams, top.idef.il,
472 fr->ePBC, fr->bMolPBC, cr, state->box);
476 if (ir->efep != efepNO)
478 /* Set free energy calculation frequency as the greatest common
479 * denominator of nstdhdl and repl_ex_nst. */
480 nstfep = ir->fepvals->nstdhdl;
483 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
485 if (useReplicaExchange)
487 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
491 /* Be REALLY careful about what flags you set here. You CANNOT assume
492 * this is the first step, since we might be restarting from a checkpoint,
493 * and in that case we should not do any modifications to the state.
495 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
497 // When restarting from a checkpoint, it can be appropriate to
498 // initialize ekind from quantities in the checkpoint. Otherwise,
499 // compute_globals must initialize ekind before the simulation
500 // starts/restarts. However, only the master rank knows what was
501 // found in the checkpoint file, so we have to communicate in
502 // order to coordinate the restart.
504 // TODO Consider removing this communication if/when checkpoint
505 // reading directly follows .tpr reading, because all ranks can
506 // agree on hasReadEkinState at that time.
507 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
510 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr);
512 if (hasReadEkinState)
514 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
517 cglo_flags = (CGLO_INITIALIZATION | CGLO_TEMPERATURE | CGLO_GSTAT
518 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
519 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
520 | (hasReadEkinState ? CGLO_READEKIN : 0));
522 bSumEkinhOld = FALSE;
524 t_vcm vcm(top_global->groups, *ir);
525 reportComRemovalInfo(fplog, vcm);
527 /* To minimize communication, compute_globals computes the COM velocity
528 * and the kinetic energy for the velocities without COM motion removed.
529 * Thus to get the kinetic energy without the COM contribution, we need
530 * to call compute_globals twice.
532 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
534 int cglo_flags_iteration = cglo_flags;
535 if (bStopCM && cgloIteration == 0)
537 cglo_flags_iteration |= CGLO_STOPCM;
538 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
540 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
541 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
542 constr, &nullSignaller, state->box,
543 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags_iteration
544 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
546 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
547 top_global, &top, state,
548 &shouldCheckNumberOfBondedInteractions);
549 if (ir->eI == eiVVAK)
551 /* a second call to get the half step temperature initialized as well */
552 /* we do the same call as above, but turn the pressure off -- internally to
553 compute_globals, this is recognized as a velocity verlet half-step
554 kinetic energy calculation. This minimized excess variables, but
555 perhaps loses some logic?*/
557 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
558 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
559 constr, &nullSignaller, state->box,
560 nullptr, &bSumEkinhOld,
561 cglo_flags & ~CGLO_PRESSURE);
564 /* Calculate the initial half step temperature, and save the ekinh_old */
565 if (startingBehavior == StartingBehavior::NewSimulation)
567 for (i = 0; (i < ir->opts.ngtc); i++)
569 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
573 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
574 temperature control */
575 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
579 if (!ir->bContinuation)
581 if (constr && ir->eConstrAlg == econtLINCS)
584 "RMS relative constraint deviation after constraining: %.2e\n",
587 if (EI_STATE_VELOCITY(ir->eI))
589 real temp = enerd->term[F_TEMP];
592 /* Result of Ekin averaged over velocities of -half
593 * and +half step, while we only have -half step here.
597 fprintf(fplog, "Initial temperature: %g K\n", temp);
602 fprintf(stderr, "starting mdrun '%s'\n",
603 *(top_global->name));
606 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
610 sprintf(tbuf, "%s", "infinite");
612 if (ir->init_step > 0)
614 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
615 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
616 gmx_step_str(ir->init_step, sbuf2),
617 ir->init_step*ir->delta_t);
621 fprintf(stderr, "%s steps, %s ps.\n",
622 gmx_step_str(ir->nsteps, sbuf), tbuf);
624 fprintf(fplog, "\n");
627 walltime_accounting_start_time(walltime_accounting);
628 wallcycle_start(wcycle, ewcRUN);
629 print_start(fplog, cr, walltime_accounting, "mdrun");
632 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
633 int chkpt_ret = fcCheckPointParallel( cr->nodeid,
637 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
641 /***********************************************************
645 ************************************************************/
648 /* Skip the first Nose-Hoover integration when we get the state from tpx */
649 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
650 bSumEkinhOld = FALSE;
652 bNeedRepartition = FALSE;
654 bool simulationsShareState = false;
655 int nstSignalComm = nstglobalcomm;
657 // TODO This implementation of ensemble orientation restraints is nasty because
658 // a user can't just do multi-sim with single-sim orientation restraints.
659 bool usingEnsembleRestraints = (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
660 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
662 // Replica exchange, ensemble restraints and AWH need all
663 // simulations to remain synchronized, so they need
664 // checkpoints and stop conditions to act on the same step, so
665 // the propagation of such signals must take place between
666 // simulations, not just within simulations.
667 // TODO: Make algorithm initializers set these flags.
668 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
670 if (simulationsShareState)
672 // Inter-simulation signal communication does not need to happen
673 // often, so we use a minimum of 200 steps to reduce overhead.
674 const int c_minimumInterSimulationSignallingInterval = 200;
675 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1)/nstglobalcomm)*nstglobalcomm;
679 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
680 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
681 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
682 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
684 auto checkpointHandler = std::make_unique<CheckpointHandler>(
685 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]),
686 simulationsShareState, ir->nstlist == 0, MASTER(cr),
687 mdrunOptions.writeConfout, mdrunOptions.checkpointOptions.period);
689 const bool resetCountersIsLocal = true;
690 auto resetHandler = std::make_unique<ResetHandler>(
691 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]), !resetCountersIsLocal,
692 ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
693 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
695 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
697 step = ir->init_step;
700 // TODO extract this to new multi-simulation module
701 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
703 if (!multisim_int_all_are_equal(ms, ir->nsteps))
705 GMX_LOG(mdlog.warning).appendText(
706 "Note: The number of steps is not consistent across multi simulations,\n"
707 "but we are proceeding anyway!");
709 if (!multisim_int_all_are_equal(ms, ir->init_step))
711 GMX_LOG(mdlog.warning).appendText(
712 "Note: The initial step is not consistent across multi simulations,\n"
713 "but we are proceeding anyway!");
717 /* and stop now if we should */
718 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
722 /* Determine if this is a neighbor search step */
723 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
725 if (bPMETune && bNStList)
727 /* PME grid + cut-off optimization with GPUs or PME nodes */
728 pme_loadbal_do(pme_loadbal, cr,
729 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
737 wallcycle_start(wcycle, ewcSTEP);
739 bLastStep = (step_rel == ir->nsteps);
740 t = t0 + step*ir->delta_t;
742 // TODO Refactor this, so that nstfep does not need a default value of zero
743 if (ir->efep != efepNO || ir->bSimTemp)
745 /* find and set the current lambdas */
746 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
748 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
749 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
750 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
751 && (ir->bExpanded) && (step > 0) &&
752 (startingBehavior == StartingBehavior::NewSimulation));
755 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep &&
756 do_per_step(step, replExParams.exchangeInterval));
758 if (doSimulatedAnnealing)
760 update_annealing_target_temp(ir, t, &upd);
763 /* Stop Center of Mass motion */
764 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
766 /* Determine whether or not to do Neighbour Searching */
767 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
769 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
771 /* do_log triggers energy and virial calculation. Because this leads
772 * to different code paths, forces can be different. Thus for exact
773 * continuation we should avoid extra log output.
774 * Note that the || bLastStep can result in non-exact continuation
775 * beyond the last step. But we don't consider that to be an issue.
778 (do_per_step(step, ir->nstlog) ||
779 (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) ||
781 do_verbose = mdrunOptions.verbose &&
782 (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
784 if (bNS && !(bFirstStep && ir->bContinuation))
786 bMasterState = FALSE;
787 /* Correct the new box if it is too skewed */
788 if (inputrecDynamicBox(ir))
790 if (correct_box(fplog, step, state->box, graph))
795 if (DOMAINDECOMP(cr) && bMasterState)
797 dd_collect_state(cr->dd, state, state_global);
800 if (DOMAINDECOMP(cr))
802 /* Repartition the domain decomposition */
803 dd_partition_system(fplog, mdlog, step, cr,
804 bMasterState, nstglobalcomm,
805 state_global, *top_global, ir, imdSession,
807 state, &f, mdAtoms, &top, fr,
810 do_verbose && !bPMETunePrinting);
811 shouldCheckNumberOfBondedInteractions = true;
812 upd.setNumAtoms(state->natoms);
816 if (MASTER(cr) && do_log)
818 energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
821 if (ir->efep != efepNO)
823 update_mdatoms(mdatoms, state->lambda[efptMASS]);
829 /* We need the kinetic energy at minus the half step for determining
830 * the full step kinetic energy and possibly for T-coupling.*/
831 /* This may not be quite working correctly yet . . . . */
832 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
833 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
834 constr, &nullSignaller, state->box,
835 &totalNumberOfBondedInteractions, &bSumEkinhOld,
836 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
837 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
838 top_global, &top, state,
839 &shouldCheckNumberOfBondedInteractions);
841 clear_mat(force_vir);
843 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
845 /* Determine the energy and pressure:
846 * at nstcalcenergy steps and at energy output steps (set below).
848 if (EI_VV(ir->eI) && (!bInitStep))
850 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
851 bCalcVir = bCalcEnerStep ||
852 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
856 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
857 bCalcVir = bCalcEnerStep ||
858 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
860 bCalcEner = bCalcEnerStep;
862 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
864 if (do_ene || do_log || bDoReplEx)
870 /* Do we need global communication ? */
871 bGStat = (bCalcVir || bCalcEner || bStopCM ||
872 do_per_step(step, nstglobalcomm) ||
873 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
875 force_flags = (GMX_FORCE_STATECHANGED |
876 ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0) |
877 GMX_FORCE_ALLFORCES |
878 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
879 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
880 (bDoFEP ? GMX_FORCE_DHDL : 0)
885 /* Now is the time to relax the shells */
886 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
887 enforcedRotation, step,
888 ir, imdSession, pull_work, bNS, force_flags, &top,
890 state, f.arrayRefWithPadding(), force_vir, mdatoms,
892 shellfc, fr, ppForceWorkload, t, mu_tot,
894 ddBalanceRegionHandler);
898 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
899 (or the AWH update will be performed twice for one step when continuing). It would be best to
900 call this update function from do_md_trajectory_writing but that would occur after do_force.
901 One would have to divide the update_awh function into one function applying the AWH force
902 and one doing the AWH bias update. The update AWH bias function could then be called after
903 do_md_trajectory_writing (then containing update_awh_history).
904 The checkpointing will in the future probably moved to the start of the md loop which will
905 rid of this issue. */
906 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
908 awh->updateHistory(state_global->awhHistory.get());
911 /* The coordinates (x) are shifted (to get whole molecules)
913 * This is parallellized as well, and does communication too.
914 * Check comments in sim_util.c
916 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation, imdSession,
918 step, nrnb, wcycle, &top,
919 state->box, state->x.arrayRefWithPadding(), &state->hist,
920 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
921 state->lambda, graph,
922 fr, ppForceWorkload, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
923 (bNS ? GMX_FORCE_NS : 0) | force_flags,
924 ddBalanceRegionHandler);
927 if (EI_VV(ir->eI) && startingBehavior == StartingBehavior::NewSimulation)
928 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
930 rvec *vbuf = nullptr;
932 wallcycle_start(wcycle, ewcUPDATE);
933 if (ir->eI == eiVV && bInitStep)
935 /* if using velocity verlet with full time step Ekin,
936 * take the first half step only to compute the
937 * virial for the first step. From there,
938 * revert back to the initial coordinates
939 * so that the input is actually the initial step.
941 snew(vbuf, state->natoms);
942 copy_rvecn(state->v.rvec_array(), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
946 /* this is for NHC in the Ekin(t+dt/2) version of vv */
947 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
950 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
951 ekind, M, &upd, etrtVELOCITY1,
954 wallcycle_stop(wcycle, ewcUPDATE);
955 constrain_velocities(step, nullptr,
959 bCalcVir, do_log, do_ene);
960 wallcycle_start(wcycle, ewcUPDATE);
961 /* if VV, compute the pressure and constraints */
962 /* For VV2, we strictly only need this if using pressure
963 * control, but we really would like to have accurate pressures
965 * Think about ways around this in the future?
966 * For now, keep this choice in comments.
968 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
969 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
971 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
972 if (bCalcEner && ir->eI == eiVVAK)
976 /* for vv, the first half of the integration actually corresponds to the previous step.
977 So we need information from the last step in the first half of the integration */
978 if (bGStat || do_per_step(step-1, nstglobalcomm))
980 wallcycle_stop(wcycle, ewcUPDATE);
981 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
982 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
983 constr, &nullSignaller, state->box,
984 &totalNumberOfBondedInteractions, &bSumEkinhOld,
985 (bGStat ? CGLO_GSTAT : 0)
986 | (bCalcEner ? CGLO_ENERGY : 0)
987 | (bTemp ? CGLO_TEMPERATURE : 0)
988 | (bPres ? CGLO_PRESSURE : 0)
989 | (bPres ? CGLO_CONSTRAINT : 0)
990 | (bStopCM ? CGLO_STOPCM : 0)
991 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
994 /* explanation of above:
995 a) We compute Ekin at the full time step
996 if 1) we are using the AveVel Ekin, and it's not the
997 initial step, or 2) if we are using AveEkin, but need the full
998 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
999 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1000 EkinAveVel because it's needed for the pressure */
1001 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1002 top_global, &top, state,
1003 &shouldCheckNumberOfBondedInteractions);
1004 wallcycle_start(wcycle, ewcUPDATE);
1006 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1011 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1012 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1014 /* TODO This is only needed when we're about to write
1015 * a checkpoint, because we use it after the restart
1016 * (in a kludge?). But what should we be doing if
1017 * the startingBehavior is NewSimulation or bInitStep are true? */
1018 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1020 copy_mat(shake_vir, state->svir_prev);
1021 copy_mat(force_vir, state->fvir_prev);
1023 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1025 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1026 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1027 enerd->term[F_EKIN] = trace(ekind->ekin);
1030 else if (bExchanged)
1032 wallcycle_stop(wcycle, ewcUPDATE);
1033 /* We need the kinetic energy at minus the half step for determining
1034 * the full step kinetic energy and possibly for T-coupling.*/
1035 /* This may not be quite working correctly yet . . . . */
1036 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
1037 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1038 constr, &nullSignaller, state->box,
1039 nullptr, &bSumEkinhOld,
1040 CGLO_GSTAT | CGLO_TEMPERATURE);
1041 wallcycle_start(wcycle, ewcUPDATE);
1044 /* if it's the initial step, we performed this first step just to get the constraint virial */
1045 if (ir->eI == eiVV && bInitStep)
1047 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1050 wallcycle_stop(wcycle, ewcUPDATE);
1053 /* compute the conserved quantity */
1056 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1059 last_ekin = enerd->term[F_EKIN];
1061 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1063 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1065 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1066 if (ir->efep != efepNO)
1068 sum_dhdl(enerd, state->lambda, ir->fepvals);
1072 /* ######## END FIRST UPDATE STEP ############## */
1073 /* ######## If doing VV, we now have v(dt) ###### */
1076 /* perform extended ensemble sampling in lambda - we don't
1077 actually move to the new state before outputting
1078 statistics, but if performing simulated tempering, we
1079 do update the velocities and the tau_t. */
1081 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, state->v.rvec_array(), mdatoms);
1082 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1085 copy_df_history(state_global->dfhist, state->dfhist);
1089 /* Now we have the energies and forces corresponding to the
1090 * coordinates at time t. We must output all of this before
1093 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1094 ir, state, state_global, observablesHistory,
1096 outf, energyOutput, ekind, f,
1097 checkpointHandler->isCheckpointingStep(),
1098 bRerunMD, bLastStep,
1099 mdrunOptions.writeConfout,
1101 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1102 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1104 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1105 if (startingBehavior != StartingBehavior::NewSimulation &&
1106 (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1108 copy_mat(state->svir_prev, shake_vir);
1109 copy_mat(state->fvir_prev, force_vir);
1112 stopHandler->setSignal();
1113 resetHandler->setSignal(walltime_accounting);
1115 if (bGStat || !PAR(cr))
1117 /* In parallel we only have to check for checkpointing in steps
1118 * where we do global communication,
1119 * otherwise the other nodes don't know.
1121 checkpointHandler->setSignal(walltime_accounting);
1124 /* ######### START SECOND UPDATE STEP ################# */
1126 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1129 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1131 gmx_bool bIfRandomize;
1132 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1133 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1134 if (constr && bIfRandomize)
1136 constrain_velocities(step, nullptr,
1140 bCalcVir, do_log, do_ene);
1143 /* Box is changed in update() when we do pressure coupling,
1144 * but we should still use the old box for energy corrections and when
1145 * writing it to the energy file, so it matches the trajectory files for
1146 * the same timestep above. Make a copy in a separate array.
1148 copy_mat(state->box, lastbox);
1152 wallcycle_start(wcycle, ewcUPDATE);
1153 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1156 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1157 /* We can only do Berendsen coupling after we have summed
1158 * the kinetic energy or virial. Since the happens
1159 * in global_state after update, we should only do it at
1160 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1165 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1166 update_pcouple_before_coordinates(fplog, step, ir, state,
1167 parrinellorahmanMu, M,
1173 /* velocity half-step update */
1174 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
1175 ekind, M, &upd, etrtVELOCITY2,
1179 /* Above, initialize just copies ekinh into ekin,
1180 * it doesn't copy position (for VV),
1181 * and entire integrator for MD.
1184 if (ir->eI == eiVVAK)
1186 /* We probably only need md->homenr, not state->natoms */
1187 if (state->natoms > cbuf_nalloc)
1189 cbuf_nalloc = state->natoms;
1190 srenew(cbuf, cbuf_nalloc);
1192 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1195 if (c_useGpuUpdateConstrain)
1197 integrator->copyCoordinatesToGpu(state->x.rvec_array());
1198 integrator->copyVelocitiesToGpu(state->v.rvec_array());
1199 integrator->copyForcesToGpu(as_rvec_array(f.data()));
1201 // This applies Leap-Frog, LINCS and SETTLE in a succession
1202 integrator->integrate(ir->delta_t, true, bCalcVir, shake_vir);
1204 integrator->copyCoordinatesFromGpu(state->x.rvec_array());
1205 integrator->copyVelocitiesFromGpu(state->v.rvec_array());
1209 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
1210 ekind, M, &upd, etrtPOSITION, cr, constr);
1212 wallcycle_stop(wcycle, ewcUPDATE);
1214 constrain_coordinates(step, &dvdl_constr, state,
1217 bCalcVir, do_log, do_ene);
1219 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state,
1220 cr, nrnb, wcycle, &upd, constr, do_log, do_ene);
1221 finish_update(ir, mdatoms,
1223 nrnb, wcycle, &upd, constr);
1226 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1228 updatePrevStepPullCom(pull_work, state);
1231 if (ir->eI == eiVVAK)
1233 /* erase F_EKIN and F_TEMP here? */
1234 /* just compute the kinetic energy at the half step to perform a trotter step */
1235 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
1236 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1237 constr, &nullSignaller, lastbox,
1238 nullptr, &bSumEkinhOld,
1239 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1241 wallcycle_start(wcycle, ewcUPDATE);
1242 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1243 /* now we know the scaling, we can compute the positions again again */
1244 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1246 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
1247 ekind, M, &upd, etrtPOSITION, cr, constr);
1248 wallcycle_stop(wcycle, ewcUPDATE);
1250 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1251 /* are the small terms in the shake_vir here due
1252 * to numerical errors, or are they important
1253 * physically? I'm thinking they are just errors, but not completely sure.
1254 * For now, will call without actually constraining, constr=NULL*/
1255 finish_update(ir, mdatoms,
1257 nrnb, wcycle, &upd, nullptr);
1261 /* this factor or 2 correction is necessary
1262 because half of the constraint force is removed
1263 in the vv step, so we have to double it. See
1264 the Redmine issue #1255. It is not yet clear
1265 if the factor of 2 is exact, or just a very
1266 good approximation, and this will be
1267 investigated. The next step is to see if this
1268 can be done adding a dhdl contribution from the
1269 rattle step, but this is somewhat more
1270 complicated with the current code. Will be
1271 investigated, hopefully for 4.6.3. However,
1272 this current solution is much better than
1273 having it completely wrong.
1275 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1279 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1282 if (vsite != nullptr)
1284 wallcycle_start(wcycle, ewcVSITECONSTR);
1285 if (graph != nullptr)
1287 shift_self(graph, state->box, state->x.rvec_array());
1289 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, state->v.rvec_array(),
1290 top.idef.iparams, top.idef.il,
1291 fr->ePBC, fr->bMolPBC, cr, state->box);
1293 if (graph != nullptr)
1295 unshift_self(graph, state->box, state->x.rvec_array());
1297 wallcycle_stop(wcycle, ewcVSITECONSTR);
1300 /* ############## IF NOT VV, Calculate globals HERE ############ */
1301 /* With Leap-Frog we can skip compute_globals at
1302 * non-communication steps, but we need to calculate
1303 * the kinetic energy one step before communication.
1306 // Organize to do inter-simulation signalling on steps if
1307 // and when algorithms require it.
1308 bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1310 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1312 // Since we're already communicating at this step, we
1313 // can propagate intra-simulation signals. Note that
1314 // check_nstglobalcomm has the responsibility for
1315 // choosing the value of nstglobalcomm that is one way
1316 // bGStat becomes true, so we can't get into a
1317 // situation where e.g. checkpointing can't be
1319 bool doIntraSimSignal = true;
1320 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1322 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, &vcm,
1323 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1326 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1327 (bGStat ? CGLO_GSTAT : 0)
1328 | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1329 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1330 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1331 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
1333 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1335 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1336 top_global, &top, state,
1337 &shouldCheckNumberOfBondedInteractions);
1341 /* ############# END CALC EKIN AND PRESSURE ################# */
1343 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1344 the virial that should probably be addressed eventually. state->veta has better properies,
1345 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1346 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1348 if (ir->efep != efepNO && !EI_VV(ir->eI))
1350 /* Sum up the foreign energy and dhdl terms for md and sd.
1351 Currently done every step so that dhdl is correct in the .edr */
1352 sum_dhdl(enerd, state->lambda, ir->fepvals);
1355 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1356 pres, force_vir, shake_vir,
1360 /* ################# END UPDATE STEP 2 ################# */
1361 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1363 /* The coordinates (x) were unshifted in update */
1366 /* We will not sum ekinh_old,
1367 * so signal that we still have to do it.
1369 bSumEkinhOld = TRUE;
1374 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1376 /* use the directly determined last velocity, not actually the averaged half steps */
1377 if (bTrotter && ir->eI == eiVV)
1379 enerd->term[F_EKIN] = last_ekin;
1381 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1383 if (integratorHasConservedEnergyQuantity(ir))
1387 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1391 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1394 /* ######### END PREPARING EDR OUTPUT ########### */
1400 if (fplog && do_log && bDoExpanded)
1402 /* only needed if doing expanded ensemble */
1403 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1404 state_global->dfhist, state->fep_state, ir->nstlog, step);
1408 energyOutput.addDataAtEnergyStep(bDoDHDL, bCalcEnerStep,
1409 t, mdatoms->tmass, enerd, state,
1410 ir->fepvals, ir->expandedvals, lastbox,
1411 shake_vir, force_vir, total_vir, pres,
1412 ekind, mu_tot, constr);
1416 energyOutput.recordNonEnergyStep();
1419 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1420 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1422 energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
1423 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1424 do_log ? fplog : nullptr,
1430 pull_print_output(pull_work, step, t);
1433 if (do_per_step(step, ir->nstlog))
1435 if (fflush(fplog) != 0)
1437 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1443 /* Have to do this part _after_ outputting the logfile and the edr file */
1444 /* Gets written into the state at the beginning of next loop*/
1445 state->fep_state = lamnew;
1447 /* Print the remaining wall clock time for the run */
1448 if (isMasterSimMasterRank(ms, MASTER(cr)) &&
1449 (do_verbose || gmx_got_usr_signal()) &&
1454 fprintf(stderr, "\n");
1456 print_time(stderr, walltime_accounting, step, ir, cr);
1459 /* Ion/water position swapping.
1460 * Not done in last step since trajectory writing happens before this call
1461 * in the MD loop and exchanges would be lost anyway. */
1462 bNeedRepartition = FALSE;
1463 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1464 do_per_step(step, ir->swap->nstswap))
1466 bNeedRepartition = do_swapcoords(cr, step, t, ir, swap, wcycle,
1467 as_rvec_array(state->x.data()),
1469 MASTER(cr) && mdrunOptions.verbose,
1472 if (bNeedRepartition && DOMAINDECOMP(cr))
1474 dd_collect_state(cr->dd, state, state_global);
1478 /* Replica exchange */
1482 bExchanged = replica_exchange(fplog, cr, ms, repl_ex,
1483 state_global, enerd,
1487 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1489 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1,
1490 state_global, *top_global, ir, imdSession,
1492 state, &f, mdAtoms, &top, fr,
1494 nrnb, wcycle, FALSE);
1495 shouldCheckNumberOfBondedInteractions = true;
1496 upd.setNumAtoms(state->natoms);
1502 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1503 /* With all integrators, except VV, we need to retain the pressure
1504 * at the current step for coupling at the next step.
1506 if ((state->flags & (1<<estPRES_PREV)) &&
1508 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1510 /* Store the pressure in t_state for pressure coupling
1511 * at the next MD step.
1513 copy_mat(pres, state->pres_prev);
1516 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1518 if ( (membed != nullptr) && (!bLastStep) )
1520 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1523 cycles = wallcycle_stop(wcycle, ewcSTEP);
1524 if (DOMAINDECOMP(cr) && wcycle)
1526 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1529 /* increase the MD step number */
1533 resetHandler->resetCounters(
1534 step, step_rel, mdlog, fplog, cr, fr->nbv.get(),
1535 nrnb, fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1537 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1538 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1541 /* End of main MD loop */
1543 /* Closing TNG files can include compressing data. Therefore it is good to do that
1544 * before stopping the time measurements. */
1545 mdoutf_tng_close(outf);
1547 /* Stop measuring walltime */
1548 walltime_accounting_end_time(walltime_accounting);
1550 if (!thisRankHasDuty(cr, DUTY_PME))
1552 /* Tell the PME only node to finish */
1553 gmx_pme_send_finish(cr);
1558 if (ir->nstcalcenergy > 0)
1560 energyOutput.printAnnealingTemperatures(fplog, groups, &(ir->opts));
1561 energyOutput.printAverages(fplog, groups);
1568 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1571 done_shellfc(fplog, shellfc, step_rel);
1573 if (useReplicaExchange && MASTER(cr))
1575 print_replica_exchange_statistics(fplog, repl_ex);
1578 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1580 global_stat_destroy(gstat);