Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / mdrun / md.cpp
blob167ad8a2fb44042350c0d0989def268467c99cd8
1 /*
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.
37 /*! \internal \file
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
44 #include "gmxpre.h"
46 #include <cinttypes>
47 #include <cmath>
48 #include <cstdio>
49 #include <cstdlib>
51 #include <algorithm>
52 #include <memory>
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 "legacysimulator.h"
141 #include "replicaexchange.h"
142 #include "shellfc.h"
144 #if GMX_FAHCORE
145 #include "corewrap.h"
146 #endif
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::LegacySimulator::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 unsigned int force_flags;
171 tensor force_vir = {{0}}, shake_vir = {{0}}, total_vir = {{0}},
172 tmp_vir = {{0}}, pres = {{0}};
173 int i, m;
174 rvec mu_tot;
175 matrix parrinellorahmanMu, M;
176 gmx_repl_ex_t repl_ex = nullptr;
177 gmx_localtop_t top;
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;
184 real dvdl_constr;
185 rvec *cbuf = nullptr;
186 int cbuf_nalloc = 0;
187 matrix lastbox;
188 int lamnew = 0;
189 /* for FEP */
190 int nstfep = 0;
191 double cycles;
192 real saved_conserved_quantity = 0;
193 real last_ekin = 0;
194 t_extmass MassQ;
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
222 // interface.
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),
245 top_global,
246 ir, cr, constr,
247 state_global, observablesHistory,
248 oenv,
249 startingBehavior);
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);
259 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, oenv, wcycle,
260 StartingBehavior::NewSimulation);
261 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work, mdoutf_get_fp_dhdl(outf), false);
263 /* Kinetic energy data */
264 std::unique_ptr<gmx_ekindata_t> eKinData = std::make_unique<gmx_ekindata_t>();
265 gmx_ekindata_t *ekind = eKinData.get();
266 init_ekindata(fplog, top_global, &(ir->opts), ekind);
267 /* Copy the cos acceleration to the groups struct */
268 ekind->cosacc.cos_accel = ir->cos_accel;
270 gstat = global_stat_init(ir);
272 /* Check for polarizable models and flexible constraints */
273 shellfc = init_shell_flexcon(fplog,
274 top_global, constr ? constr->numFlexibleConstraints() : 0,
275 ir->nstcalcenergy, DOMAINDECOMP(cr));
278 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
279 if ((io > 2000) && MASTER(cr))
281 fprintf(stderr,
282 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
283 io);
287 // Local state only becomes valid now.
288 std::unique_ptr<t_state> stateInstance;
289 t_state * state;
292 auto mdatoms = mdAtoms->mdatoms();
294 std::unique_ptr<UpdateConstrainCuda> integrator;
296 if (DOMAINDECOMP(cr))
298 GMX_RELEASE_ASSERT(!c_useGpuUpdateConstrain, "Domain decomposition is not supported with GPU-based update-constraints.");
299 dd_init_local_top(*top_global, &top);
301 stateInstance = std::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, imdSession,
308 pull_work,
309 state, &f, mdAtoms, &top, fr,
310 vsite, constr,
311 nrnb, nullptr, FALSE);
312 shouldCheckNumberOfBondedInteractions = true;
313 upd.setNumAtoms(state->natoms);
315 else
317 state_change_natoms(state_global, state_global->natoms);
318 f.resizeWithPadding(state_global->natoms);
319 /* Copy the pointer to the global state */
320 state = state_global;
322 /* Generate and initialize new topology */
323 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr,
324 &graph, mdAtoms, constr, vsite, shellfc);
326 upd.setNumAtoms(state->natoms);
328 if (c_useGpuUpdateConstrain)
330 GMX_RELEASE_ASSERT(ir->eI == eiMD, "Only md integrator is supported on the GPU.");
331 GMX_RELEASE_ASSERT(ir->etc == etcNO, "Temperature coupling is not supported on the GPU.");
332 GMX_RELEASE_ASSERT(ir->epc == epcNO, "Pressure coupling is not supported on the GPU.");
333 GMX_RELEASE_ASSERT(!mdatoms->haveVsites, "Virtual sites are not supported on the GPU");
334 GMX_LOG(mdlog.info).asParagraph().
335 appendText("Using CUDA GPU-based update and constraints module.");
336 integrator = std::make_unique<UpdateConstrainCuda>(state->natoms, *ir, *top_global);
337 integrator->set(top.idef, *mdatoms);
338 t_pbc pbc;
339 set_pbc(&pbc, epbcXYZ, state->box);
340 integrator->setPbc(&pbc);
345 if (fr->nbv->useGpu())
347 changePinningPolicy(&state->x, gmx::PinningPolicy::PinnedIfSupported);
350 // NOTE: The global state is no longer used at this point.
351 // But state_global is still used as temporary storage space for writing
352 // the global state to file and potentially for replica exchange.
353 // (Global topology should persist.)
355 update_mdatoms(mdatoms, state->lambda[efptMASS]);
357 if (ir->bExpanded)
359 /* Check nstexpanded here, because the grompp check was broken */
360 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
362 gmx_fatal(FARGS, "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
364 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation,
365 ir, state->dfhist);
368 if (MASTER(cr))
370 if (startingBehavior != StartingBehavior::NewSimulation)
372 /* Restore from energy history if appending to output files */
373 if (startingBehavior == StartingBehavior::RestartWithAppending)
375 /* If no history is available (because a checkpoint is from before
376 * it was written) make a new one later, otherwise restore it.
378 if (observablesHistory->energyHistory)
380 energyOutput.restoreFromEnergyHistory(*observablesHistory->energyHistory);
383 else if (observablesHistory->energyHistory)
385 /* We might have read an energy history from checkpoint.
386 * As we are not appending, we want to restart the statistics.
387 * Free the allocated memory and reset the counts.
389 observablesHistory->energyHistory = {};
390 /* We might have read a pull history from checkpoint.
391 * We will still want to keep the statistics, so that the files
392 * can be joined and still be meaningful.
393 * This means that observablesHistory->pullHistory
394 * should not be reset.
398 if (!observablesHistory->energyHistory)
400 observablesHistory->energyHistory = std::make_unique<energyhistory_t>();
402 if (!observablesHistory->pullHistory)
404 observablesHistory->pullHistory = std::make_unique<PullHistory>();
406 /* Set the initial energy history */
407 energyOutput.fillEnergyHistory(observablesHistory->energyHistory.get());
410 preparePrevStepPullCom(ir, pull_work, mdatoms, state, state_global, cr,
411 startingBehavior != StartingBehavior::NewSimulation);
413 // TODO: Remove this by converting AWH into a ForceProvider
414 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms, startingBehavior != StartingBehavior::NewSimulation,
415 shellfc != nullptr,
416 opt2fn("-awh", nfile, fnm), pull_work);
418 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
419 if (useReplicaExchange && MASTER(cr))
421 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir,
422 replExParams);
424 /* PME tuning is only supported in the Verlet scheme, with PME for
425 * Coulomb. It is not supported with only LJ PME. */
426 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) &&
427 !mdrunOptions.reproducible && ir->cutoff_scheme != ecutsGROUP);
429 pme_load_balancing_t *pme_loadbal = nullptr;
430 if (bPMETune)
432 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box,
433 *fr->ic, *fr->nbv, fr->pmedata, fr->nbv->useGpu(),
434 &bPMETunePrinting);
437 if (!ir->bContinuation)
439 if (state->flags & (1u << estV))
441 auto v = makeArrayRef(state->v);
442 /* Set the velocities of vsites, shells and frozen atoms to zero */
443 for (i = 0; i < mdatoms->homenr; i++)
445 if (mdatoms->ptype[i] == eptVSite ||
446 mdatoms->ptype[i] == eptShell)
448 clear_rvec(v[i]);
450 else if (mdatoms->cFREEZE)
452 for (m = 0; m < DIM; m++)
454 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
456 v[i][m] = 0;
463 if (constr)
465 /* Constrain the initial coordinates and velocities */
466 do_constrain_first(fplog, constr, ir, mdatoms, state);
468 if (vsite)
470 /* Construct the virtual sites for the initial configuration */
471 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, nullptr,
472 top.idef.iparams, top.idef.il,
473 fr->ePBC, fr->bMolPBC, cr, state->box);
477 if (ir->efep != efepNO)
479 /* Set free energy calculation frequency as the greatest common
480 * denominator of nstdhdl and repl_ex_nst. */
481 nstfep = ir->fepvals->nstdhdl;
482 if (ir->bExpanded)
484 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
486 if (useReplicaExchange)
488 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
492 /* Be REALLY careful about what flags you set here. You CANNOT assume
493 * this is the first step, since we might be restarting from a checkpoint,
494 * and in that case we should not do any modifications to the state.
496 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
498 // When restarting from a checkpoint, it can be appropriate to
499 // initialize ekind from quantities in the checkpoint. Otherwise,
500 // compute_globals must initialize ekind before the simulation
501 // starts/restarts. However, only the master rank knows what was
502 // found in the checkpoint file, so we have to communicate in
503 // order to coordinate the restart.
505 // TODO Consider removing this communication if/when checkpoint
506 // reading directly follows .tpr reading, because all ranks can
507 // agree on hasReadEkinState at that time.
508 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
509 if (PAR(cr))
511 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr);
513 if (hasReadEkinState)
515 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
518 unsigned int cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
519 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
520 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
521 | (hasReadEkinState ? CGLO_READEKIN : 0));
523 bSumEkinhOld = FALSE;
525 t_vcm vcm(top_global->groups, *ir);
526 reportComRemovalInfo(fplog, vcm);
528 /* To minimize communication, compute_globals computes the COM velocity
529 * and the kinetic energy for the velocities without COM motion removed.
530 * Thus to get the kinetic energy without the COM contribution, we need
531 * to call compute_globals twice.
533 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
535 unsigned int cglo_flags_iteration = cglo_flags;
536 if (bStopCM && cgloIteration == 0)
538 cglo_flags_iteration |= CGLO_STOPCM;
539 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
541 compute_globals(gstat, cr, ir, fr, ekind,
542 state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
543 mdatoms, nrnb, &vcm,
544 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
545 constr, &nullSignaller, state->box,
546 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags_iteration
547 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
548 if (cglo_flags_iteration & CGLO_STOPCM)
550 /* At initialization, do not pass x with acceleration-correction mode
551 * to avoid (incorrect) correction of the initial coordinates.
553 rvec *xPtr = nullptr;
554 if (vcm.mode != ecmLINEAR_ACCELERATION_CORRECTION)
556 xPtr = state->x.rvec_array();
558 process_and_stopcm_grp(fplog, &vcm, *mdatoms, xPtr, state->v.rvec_array());
559 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
562 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
563 top_global, &top, state->x.rvec_array(), state->box,
564 &shouldCheckNumberOfBondedInteractions);
565 if (ir->eI == eiVVAK)
567 /* a second call to get the half step temperature initialized as well */
568 /* we do the same call as above, but turn the pressure off -- internally to
569 compute_globals, this is recognized as a velocity verlet half-step
570 kinetic energy calculation. This minimized excess variables, but
571 perhaps loses some logic?*/
573 compute_globals(gstat, cr, ir, fr, ekind,
574 state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
575 mdatoms, nrnb, &vcm,
576 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
577 constr, &nullSignaller, state->box,
578 nullptr, &bSumEkinhOld,
579 cglo_flags & ~CGLO_PRESSURE);
582 /* Calculate the initial half step temperature, and save the ekinh_old */
583 if (startingBehavior == StartingBehavior::NewSimulation)
585 for (i = 0; (i < ir->opts.ngtc); i++)
587 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
591 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
592 temperature control */
593 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
595 if (MASTER(cr))
597 if (!ir->bContinuation)
599 if (constr && ir->eConstrAlg == econtLINCS)
601 fprintf(fplog,
602 "RMS relative constraint deviation after constraining: %.2e\n",
603 constr->rmsd());
605 if (EI_STATE_VELOCITY(ir->eI))
607 real temp = enerd->term[F_TEMP];
608 if (ir->eI != eiVV)
610 /* Result of Ekin averaged over velocities of -half
611 * and +half step, while we only have -half step here.
613 temp *= 2;
615 fprintf(fplog, "Initial temperature: %g K\n", temp);
619 char tbuf[20];
620 fprintf(stderr, "starting mdrun '%s'\n",
621 *(top_global->name));
622 if (ir->nsteps >= 0)
624 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
626 else
628 sprintf(tbuf, "%s", "infinite");
630 if (ir->init_step > 0)
632 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
633 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
634 gmx_step_str(ir->init_step, sbuf2),
635 ir->init_step*ir->delta_t);
637 else
639 fprintf(stderr, "%s steps, %s ps.\n",
640 gmx_step_str(ir->nsteps, sbuf), tbuf);
642 fprintf(fplog, "\n");
645 walltime_accounting_start_time(walltime_accounting);
646 wallcycle_start(wcycle, ewcRUN);
647 print_start(fplog, cr, walltime_accounting, "mdrun");
649 #if GMX_FAHCORE
650 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
651 int chkpt_ret = fcCheckPointParallel( cr->nodeid,
652 NULL, 0);
653 if (chkpt_ret == 0)
655 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
657 #endif
659 /***********************************************************
661 * Loop over MD steps
663 ************************************************************/
665 bFirstStep = TRUE;
666 /* Skip the first Nose-Hoover integration when we get the state from tpx */
667 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
668 bSumEkinhOld = FALSE;
669 bExchanged = FALSE;
670 bNeedRepartition = FALSE;
672 bool simulationsShareState = false;
673 int nstSignalComm = nstglobalcomm;
675 // TODO This implementation of ensemble orientation restraints is nasty because
676 // a user can't just do multi-sim with single-sim orientation restraints.
677 bool usingEnsembleRestraints = (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
678 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
680 // Replica exchange, ensemble restraints and AWH need all
681 // simulations to remain synchronized, so they need
682 // checkpoints and stop conditions to act on the same step, so
683 // the propagation of such signals must take place between
684 // simulations, not just within simulations.
685 // TODO: Make algorithm initializers set these flags.
686 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
688 if (simulationsShareState)
690 // Inter-simulation signal communication does not need to happen
691 // often, so we use a minimum of 200 steps to reduce overhead.
692 const int c_minimumInterSimulationSignallingInterval = 200;
693 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1)/nstglobalcomm)*nstglobalcomm;
697 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
698 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
699 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
700 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
702 auto checkpointHandler = std::make_unique<CheckpointHandler>(
703 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]),
704 simulationsShareState, ir->nstlist == 0, MASTER(cr),
705 mdrunOptions.writeConfout, mdrunOptions.checkpointOptions.period);
707 const bool resetCountersIsLocal = true;
708 auto resetHandler = std::make_unique<ResetHandler>(
709 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]), !resetCountersIsLocal,
710 ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
711 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
713 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
715 step = ir->init_step;
716 step_rel = 0;
718 // TODO extract this to new multi-simulation module
719 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
721 if (!multisim_int_all_are_equal(ms, ir->nsteps))
723 GMX_LOG(mdlog.warning).appendText(
724 "Note: The number of steps is not consistent across multi simulations,\n"
725 "but we are proceeding anyway!");
727 if (!multisim_int_all_are_equal(ms, ir->init_step))
729 GMX_LOG(mdlog.warning).appendText(
730 "Note: The initial step is not consistent across multi simulations,\n"
731 "but we are proceeding anyway!");
735 /* and stop now if we should */
736 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
737 while (!bLastStep)
740 /* Determine if this is a neighbor search step */
741 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
743 if (bPMETune && bNStList)
745 /* PME grid + cut-off optimization with GPUs or PME nodes */
746 pme_loadbal_do(pme_loadbal, cr,
747 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
748 fplog, mdlog,
749 *ir, fr, *state,
750 wcycle,
751 step, step_rel,
752 &bPMETunePrinting);
755 wallcycle_start(wcycle, ewcSTEP);
757 bLastStep = (step_rel == ir->nsteps);
758 t = t0 + step*ir->delta_t;
760 // TODO Refactor this, so that nstfep does not need a default value of zero
761 if (ir->efep != efepNO || ir->bSimTemp)
763 /* find and set the current lambdas */
764 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
766 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
767 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
768 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
769 && (ir->bExpanded) && (step > 0) &&
770 (startingBehavior == StartingBehavior::NewSimulation));
773 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep &&
774 do_per_step(step, replExParams.exchangeInterval));
776 if (doSimulatedAnnealing)
778 update_annealing_target_temp(ir, t, &upd);
781 /* Stop Center of Mass motion */
782 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
784 /* Determine whether or not to do Neighbour Searching */
785 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
787 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
789 /* do_log triggers energy and virial calculation. Because this leads
790 * to different code paths, forces can be different. Thus for exact
791 * continuation we should avoid extra log output.
792 * Note that the || bLastStep can result in non-exact continuation
793 * beyond the last step. But we don't consider that to be an issue.
795 do_log =
796 (do_per_step(step, ir->nstlog) ||
797 (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) ||
798 bLastStep);
799 do_verbose = mdrunOptions.verbose &&
800 (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
802 if (bNS && !(bFirstStep && ir->bContinuation))
804 bMasterState = FALSE;
805 /* Correct the new box if it is too skewed */
806 if (inputrecDynamicBox(ir))
808 if (correct_box(fplog, step, state->box, graph))
810 bMasterState = TRUE;
813 if (DOMAINDECOMP(cr) && bMasterState)
815 dd_collect_state(cr->dd, state, state_global);
818 if (DOMAINDECOMP(cr))
820 /* Repartition the domain decomposition */
821 dd_partition_system(fplog, mdlog, step, cr,
822 bMasterState, nstglobalcomm,
823 state_global, *top_global, ir, imdSession,
824 pull_work,
825 state, &f, mdAtoms, &top, fr,
826 vsite, constr,
827 nrnb, wcycle,
828 do_verbose && !bPMETunePrinting);
829 shouldCheckNumberOfBondedInteractions = true;
830 upd.setNumAtoms(state->natoms);
834 if (MASTER(cr) && do_log)
836 energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
839 if (ir->efep != efepNO)
841 update_mdatoms(mdatoms, state->lambda[efptMASS]);
844 if (bExchanged)
847 /* We need the kinetic energy at minus the half step for determining
848 * the full step kinetic energy and possibly for T-coupling.*/
849 /* This may not be quite working correctly yet . . . . */
850 compute_globals(gstat, cr, ir, fr, ekind,
851 state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
852 mdatoms, nrnb, &vcm,
853 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
854 constr, &nullSignaller, state->box,
855 &totalNumberOfBondedInteractions, &bSumEkinhOld,
856 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
857 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
858 top_global, &top, state->x.rvec_array(), state->box,
859 &shouldCheckNumberOfBondedInteractions);
861 clear_mat(force_vir);
863 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
865 /* Determine the energy and pressure:
866 * at nstcalcenergy steps and at energy output steps (set below).
868 if (EI_VV(ir->eI) && (!bInitStep))
870 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
871 bCalcVir = bCalcEnerStep ||
872 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
874 else
876 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
877 bCalcVir = bCalcEnerStep ||
878 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
880 bCalcEner = bCalcEnerStep;
882 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
884 if (do_ene || do_log || bDoReplEx)
886 bCalcVir = TRUE;
887 bCalcEner = TRUE;
890 /* Do we need global communication ? */
891 bGStat = (bCalcVir || bCalcEner || bStopCM ||
892 do_per_step(step, nstglobalcomm) ||
893 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
895 force_flags = (GMX_FORCE_STATECHANGED |
896 ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0) |
897 GMX_FORCE_ALLFORCES |
898 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
899 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
900 (bDoFEP ? GMX_FORCE_DHDL : 0)
903 if (shellfc)
905 /* Now is the time to relax the shells */
906 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
907 enforcedRotation, step,
908 ir, imdSession, pull_work, bNS, force_flags, &top,
909 constr, enerd, fcd,
910 state->natoms,
911 state->x.arrayRefWithPadding(),
912 state->v.arrayRefWithPadding(),
913 state->box,
914 state->lambda,
915 &state->hist,
916 f.arrayRefWithPadding(), force_vir, mdatoms,
917 nrnb, wcycle, graph,
918 shellfc, fr, ppForceWorkload, t, mu_tot,
919 vsite,
920 ddBalanceRegionHandler);
922 else
924 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
925 (or the AWH update will be performed twice for one step when continuing). It would be best to
926 call this update function from do_md_trajectory_writing but that would occur after do_force.
927 One would have to divide the update_awh function into one function applying the AWH force
928 and one doing the AWH bias update. The update AWH bias function could then be called after
929 do_md_trajectory_writing (then containing update_awh_history).
930 The checkpointing will in the future probably moved to the start of the md loop which will
931 rid of this issue. */
932 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
934 awh->updateHistory(state_global->awhHistory.get());
937 /* The coordinates (x) are shifted (to get whole molecules)
938 * in do_force.
939 * This is parallellized as well, and does communication too.
940 * Check comments in sim_util.c
942 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation, imdSession,
943 pull_work,
944 step, nrnb, wcycle, &top,
945 state->box, state->x.arrayRefWithPadding(), &state->hist,
946 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
947 state->lambda, graph,
948 fr, ppForceWorkload, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
949 (bNS ? GMX_FORCE_NS : 0) | force_flags,
950 ddBalanceRegionHandler);
953 if (EI_VV(ir->eI) && startingBehavior == StartingBehavior::NewSimulation)
954 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
956 rvec *vbuf = nullptr;
958 wallcycle_start(wcycle, ewcUPDATE);
959 if (ir->eI == eiVV && bInitStep)
961 /* if using velocity verlet with full time step Ekin,
962 * take the first half step only to compute the
963 * virial for the first step. From there,
964 * revert back to the initial coordinates
965 * so that the input is actually the initial step.
967 snew(vbuf, state->natoms);
968 copy_rvecn(state->v.rvec_array(), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
970 else
972 /* this is for NHC in the Ekin(t+dt/2) version of vv */
973 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
976 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
977 ekind, M, &upd, etrtVELOCITY1,
978 cr, constr);
980 wallcycle_stop(wcycle, ewcUPDATE);
981 constrain_velocities(step, nullptr,
982 state,
983 shake_vir,
984 constr,
985 bCalcVir, do_log, do_ene);
986 wallcycle_start(wcycle, ewcUPDATE);
987 /* if VV, compute the pressure and constraints */
988 /* For VV2, we strictly only need this if using pressure
989 * control, but we really would like to have accurate pressures
990 * printed out.
991 * Think about ways around this in the future?
992 * For now, keep this choice in comments.
994 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
995 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
996 bPres = TRUE;
997 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
998 if (bCalcEner && ir->eI == eiVVAK)
1000 bSumEkinhOld = TRUE;
1002 /* for vv, the first half of the integration actually corresponds to the previous step.
1003 So we need information from the last step in the first half of the integration */
1004 if (bGStat || do_per_step(step-1, nstglobalcomm))
1006 wallcycle_stop(wcycle, ewcUPDATE);
1007 compute_globals(gstat, cr, ir, fr, ekind,
1008 state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
1009 mdatoms, nrnb, &vcm,
1010 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1011 constr, &nullSignaller, state->box,
1012 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1013 (bGStat ? CGLO_GSTAT : 0)
1014 | (bCalcEner ? CGLO_ENERGY : 0)
1015 | (bTemp ? CGLO_TEMPERATURE : 0)
1016 | (bPres ? CGLO_PRESSURE : 0)
1017 | (bPres ? CGLO_CONSTRAINT : 0)
1018 | (bStopCM ? CGLO_STOPCM : 0)
1019 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1020 | CGLO_SCALEEKIN
1022 /* explanation of above:
1023 a) We compute Ekin at the full time step
1024 if 1) we are using the AveVel Ekin, and it's not the
1025 initial step, or 2) if we are using AveEkin, but need the full
1026 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1027 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1028 EkinAveVel because it's needed for the pressure */
1029 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1030 top_global, &top, state->x.rvec_array(), state->box,
1031 &shouldCheckNumberOfBondedInteractions);
1032 if (bStopCM)
1034 process_and_stopcm_grp(fplog, &vcm, *mdatoms, state->x.rvec_array(), state->v.rvec_array());
1035 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1037 wallcycle_start(wcycle, ewcUPDATE);
1039 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1040 if (!bInitStep)
1042 if (bTrotter)
1044 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1045 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1047 /* TODO This is only needed when we're about to write
1048 * a checkpoint, because we use it after the restart
1049 * (in a kludge?). But what should we be doing if
1050 * the startingBehavior is NewSimulation or bInitStep are true? */
1051 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1053 copy_mat(shake_vir, state->svir_prev);
1054 copy_mat(force_vir, state->fvir_prev);
1056 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1058 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1059 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1060 enerd->term[F_EKIN] = trace(ekind->ekin);
1063 else if (bExchanged)
1065 wallcycle_stop(wcycle, ewcUPDATE);
1066 /* We need the kinetic energy at minus the half step for determining
1067 * the full step kinetic energy and possibly for T-coupling.*/
1068 /* This may not be quite working correctly yet . . . . */
1069 compute_globals(gstat, cr, ir, fr, ekind,
1070 state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
1071 mdatoms, nrnb, &vcm,
1072 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1073 constr, &nullSignaller, state->box,
1074 nullptr, &bSumEkinhOld,
1075 CGLO_GSTAT | CGLO_TEMPERATURE);
1076 wallcycle_start(wcycle, ewcUPDATE);
1079 /* if it's the initial step, we performed this first step just to get the constraint virial */
1080 if (ir->eI == eiVV && bInitStep)
1082 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1083 sfree(vbuf);
1085 wallcycle_stop(wcycle, ewcUPDATE);
1088 /* compute the conserved quantity */
1089 if (EI_VV(ir->eI))
1091 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1092 if (ir->eI == eiVV)
1094 last_ekin = enerd->term[F_EKIN];
1096 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1098 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1100 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1101 if (ir->efep != efepNO)
1103 sum_dhdl(enerd, state->lambda, ir->fepvals);
1107 /* ######## END FIRST UPDATE STEP ############## */
1108 /* ######## If doing VV, we now have v(dt) ###### */
1109 if (bDoExpanded)
1111 /* perform extended ensemble sampling in lambda - we don't
1112 actually move to the new state before outputting
1113 statistics, but if performing simulated tempering, we
1114 do update the velocities and the tau_t. */
1116 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, state->v.rvec_array(), mdatoms);
1117 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1118 if (MASTER(cr))
1120 copy_df_history(state_global->dfhist, state->dfhist);
1124 /* Now we have the energies and forces corresponding to the
1125 * coordinates at time t. We must output all of this before
1126 * the update.
1128 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1129 ir, state, state_global, observablesHistory,
1130 top_global, fr,
1131 outf, energyOutput, ekind, f,
1132 checkpointHandler->isCheckpointingStep(),
1133 bRerunMD, bLastStep,
1134 mdrunOptions.writeConfout,
1135 bSumEkinhOld);
1136 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1137 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1139 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1140 if (startingBehavior != StartingBehavior::NewSimulation &&
1141 (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1143 copy_mat(state->svir_prev, shake_vir);
1144 copy_mat(state->fvir_prev, force_vir);
1147 stopHandler->setSignal();
1148 resetHandler->setSignal(walltime_accounting);
1150 if (bGStat || !PAR(cr))
1152 /* In parallel we only have to check for checkpointing in steps
1153 * where we do global communication,
1154 * otherwise the other nodes don't know.
1156 checkpointHandler->setSignal(walltime_accounting);
1159 /* ######### START SECOND UPDATE STEP ################# */
1161 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1162 in preprocessing */
1164 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1166 gmx_bool bIfRandomize;
1167 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1168 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1169 if (constr && bIfRandomize)
1171 constrain_velocities(step, nullptr,
1172 state,
1173 tmp_vir,
1174 constr,
1175 bCalcVir, do_log, do_ene);
1178 /* Box is changed in update() when we do pressure coupling,
1179 * but we should still use the old box for energy corrections and when
1180 * writing it to the energy file, so it matches the trajectory files for
1181 * the same timestep above. Make a copy in a separate array.
1183 copy_mat(state->box, lastbox);
1185 dvdl_constr = 0;
1187 wallcycle_start(wcycle, ewcUPDATE);
1188 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1189 if (bTrotter)
1191 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1192 /* We can only do Berendsen coupling after we have summed
1193 * the kinetic energy or virial. Since the happens
1194 * in global_state after update, we should only do it at
1195 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1198 else
1200 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1201 update_pcouple_before_coordinates(fplog, step, ir, state,
1202 parrinellorahmanMu, M,
1203 bInitStep);
1206 if (EI_VV(ir->eI))
1208 /* velocity half-step update */
1209 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
1210 ekind, M, &upd, etrtVELOCITY2,
1211 cr, constr);
1214 /* Above, initialize just copies ekinh into ekin,
1215 * it doesn't copy position (for VV),
1216 * and entire integrator for MD.
1219 if (ir->eI == eiVVAK)
1221 /* We probably only need md->homenr, not state->natoms */
1222 if (state->natoms > cbuf_nalloc)
1224 cbuf_nalloc = state->natoms;
1225 srenew(cbuf, cbuf_nalloc);
1227 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1230 if (c_useGpuUpdateConstrain)
1232 integrator->copyCoordinatesToGpu(state->x.rvec_array());
1233 integrator->copyVelocitiesToGpu(state->v.rvec_array());
1234 integrator->copyForcesToGpu(as_rvec_array(f.data()));
1236 // This applies Leap-Frog, LINCS and SETTLE in a succession
1237 integrator->integrate(ir->delta_t, true, bCalcVir, shake_vir);
1239 integrator->copyCoordinatesFromGpu(state->x.rvec_array());
1240 integrator->copyVelocitiesFromGpu(state->v.rvec_array());
1242 else
1244 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
1245 ekind, M, &upd, etrtPOSITION, cr, constr);
1247 wallcycle_stop(wcycle, ewcUPDATE);
1249 constrain_coordinates(step, &dvdl_constr, state,
1250 shake_vir,
1251 &upd, constr,
1252 bCalcVir, do_log, do_ene);
1254 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state,
1255 cr, nrnb, wcycle, &upd, constr, do_log, do_ene);
1256 finish_update(ir, mdatoms,
1257 state, graph,
1258 nrnb, wcycle, &upd, constr);
1261 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1263 updatePrevStepPullCom(pull_work, state);
1266 if (ir->eI == eiVVAK)
1268 /* erase F_EKIN and F_TEMP here? */
1269 /* just compute the kinetic energy at the half step to perform a trotter step */
1270 compute_globals(gstat, cr, ir, fr, ekind,
1271 state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
1272 mdatoms, nrnb, &vcm,
1273 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1274 constr, &nullSignaller, lastbox,
1275 nullptr, &bSumEkinhOld,
1276 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1278 wallcycle_start(wcycle, ewcUPDATE);
1279 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1280 /* now we know the scaling, we can compute the positions again again */
1281 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1283 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd,
1284 ekind, M, &upd, etrtPOSITION, cr, constr);
1285 wallcycle_stop(wcycle, ewcUPDATE);
1287 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1288 /* are the small terms in the shake_vir here due
1289 * to numerical errors, or are they important
1290 * physically? I'm thinking they are just errors, but not completely sure.
1291 * For now, will call without actually constraining, constr=NULL*/
1292 finish_update(ir, mdatoms,
1293 state, graph,
1294 nrnb, wcycle, &upd, nullptr);
1296 if (EI_VV(ir->eI))
1298 /* this factor or 2 correction is necessary
1299 because half of the constraint force is removed
1300 in the vv step, so we have to double it. See
1301 the Redmine issue #1255. It is not yet clear
1302 if the factor of 2 is exact, or just a very
1303 good approximation, and this will be
1304 investigated. The next step is to see if this
1305 can be done adding a dhdl contribution from the
1306 rattle step, but this is somewhat more
1307 complicated with the current code. Will be
1308 investigated, hopefully for 4.6.3. However,
1309 this current solution is much better than
1310 having it completely wrong.
1312 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1314 else
1316 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1319 if (vsite != nullptr)
1321 wallcycle_start(wcycle, ewcVSITECONSTR);
1322 if (graph != nullptr)
1324 shift_self(graph, state->box, state->x.rvec_array());
1326 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, state->v.rvec_array(),
1327 top.idef.iparams, top.idef.il,
1328 fr->ePBC, fr->bMolPBC, cr, state->box);
1330 if (graph != nullptr)
1332 unshift_self(graph, state->box, state->x.rvec_array());
1334 wallcycle_stop(wcycle, ewcVSITECONSTR);
1337 /* ############## IF NOT VV, Calculate globals HERE ############ */
1338 /* With Leap-Frog we can skip compute_globals at
1339 * non-communication steps, but we need to calculate
1340 * the kinetic energy one step before communication.
1343 // Organize to do inter-simulation signalling on steps if
1344 // and when algorithms require it.
1345 bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1347 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1349 // Since we're already communicating at this step, we
1350 // can propagate intra-simulation signals. Note that
1351 // check_nstglobalcomm has the responsibility for
1352 // choosing the value of nstglobalcomm that is one way
1353 // bGStat becomes true, so we can't get into a
1354 // situation where e.g. checkpointing can't be
1355 // signalled.
1356 bool doIntraSimSignal = true;
1357 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1359 compute_globals(gstat, cr, ir, fr, ekind,
1360 state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
1361 mdatoms, nrnb, &vcm,
1362 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1363 constr, &signaller,
1364 lastbox,
1365 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1366 (bGStat ? CGLO_GSTAT : 0)
1367 | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1368 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1369 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1370 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
1371 | CGLO_CONSTRAINT
1372 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1374 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1375 top_global, &top, state->x.rvec_array(), state->box,
1376 &shouldCheckNumberOfBondedInteractions);
1377 if (!EI_VV(ir->eI) && bStopCM)
1379 process_and_stopcm_grp(fplog, &vcm, *mdatoms, state->x.rvec_array(), state->v.rvec_array());
1380 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1385 /* ############# END CALC EKIN AND PRESSURE ################# */
1387 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1388 the virial that should probably be addressed eventually. state->veta has better properies,
1389 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1390 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1392 if (ir->efep != efepNO && !EI_VV(ir->eI))
1394 /* Sum up the foreign energy and dhdl terms for md and sd.
1395 Currently done every step so that dhdl is correct in the .edr */
1396 sum_dhdl(enerd, state->lambda, ir->fepvals);
1399 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1400 pres, force_vir, shake_vir,
1401 parrinellorahmanMu,
1402 state, nrnb, &upd);
1404 /* ################# END UPDATE STEP 2 ################# */
1405 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1407 /* The coordinates (x) were unshifted in update */
1408 if (!bGStat)
1410 /* We will not sum ekinh_old,
1411 * so signal that we still have to do it.
1413 bSumEkinhOld = TRUE;
1416 if (bCalcEner)
1418 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1420 /* use the directly determined last velocity, not actually the averaged half steps */
1421 if (bTrotter && ir->eI == eiVV)
1423 enerd->term[F_EKIN] = last_ekin;
1425 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1427 if (integratorHasConservedEnergyQuantity(ir))
1429 if (EI_VV(ir->eI))
1431 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1433 else
1435 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1438 /* ######### END PREPARING EDR OUTPUT ########### */
1441 /* Output stuff */
1442 if (MASTER(cr))
1444 if (fplog && do_log && bDoExpanded)
1446 /* only needed if doing expanded ensemble */
1447 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1448 state_global->dfhist, state->fep_state, ir->nstlog, step);
1450 if (bCalcEner)
1452 energyOutput.addDataAtEnergyStep(bDoDHDL, bCalcEnerStep,
1453 t, mdatoms->tmass, enerd, state,
1454 ir->fepvals, ir->expandedvals, lastbox,
1455 shake_vir, force_vir, total_vir, pres,
1456 ekind, mu_tot, constr);
1458 else
1460 energyOutput.recordNonEnergyStep();
1463 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1464 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1466 energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
1467 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1468 do_log ? fplog : nullptr,
1469 step, t,
1470 fcd, awh.get());
1472 if (ir->bPull)
1474 pull_print_output(pull_work, step, t);
1477 if (do_per_step(step, ir->nstlog))
1479 if (fflush(fplog) != 0)
1481 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1485 if (bDoExpanded)
1487 /* Have to do this part _after_ outputting the logfile and the edr file */
1488 /* Gets written into the state at the beginning of next loop*/
1489 state->fep_state = lamnew;
1491 /* Print the remaining wall clock time for the run */
1492 if (isMasterSimMasterRank(ms, MASTER(cr)) &&
1493 (do_verbose || gmx_got_usr_signal()) &&
1494 !bPMETunePrinting)
1496 if (shellfc)
1498 fprintf(stderr, "\n");
1500 print_time(stderr, walltime_accounting, step, ir, cr);
1503 /* Ion/water position swapping.
1504 * Not done in last step since trajectory writing happens before this call
1505 * in the MD loop and exchanges would be lost anyway. */
1506 bNeedRepartition = FALSE;
1507 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1508 do_per_step(step, ir->swap->nstswap))
1510 bNeedRepartition = do_swapcoords(cr, step, t, ir, swap, wcycle,
1511 as_rvec_array(state->x.data()),
1512 state->box,
1513 MASTER(cr) && mdrunOptions.verbose,
1514 bRerunMD);
1516 if (bNeedRepartition && DOMAINDECOMP(cr))
1518 dd_collect_state(cr->dd, state, state_global);
1522 /* Replica exchange */
1523 bExchanged = FALSE;
1524 if (bDoReplEx)
1526 bExchanged = replica_exchange(fplog, cr, ms, repl_ex,
1527 state_global, enerd,
1528 state, step, t);
1531 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1533 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1,
1534 state_global, *top_global, ir, imdSession,
1535 pull_work,
1536 state, &f, mdAtoms, &top, fr,
1537 vsite, constr,
1538 nrnb, wcycle, FALSE);
1539 shouldCheckNumberOfBondedInteractions = true;
1540 upd.setNumAtoms(state->natoms);
1543 bFirstStep = FALSE;
1544 bInitStep = FALSE;
1546 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1547 /* With all integrators, except VV, we need to retain the pressure
1548 * at the current step for coupling at the next step.
1550 if ((state->flags & (1u<<estPRES_PREV)) &&
1551 (bGStatEveryStep ||
1552 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1554 /* Store the pressure in t_state for pressure coupling
1555 * at the next MD step.
1557 copy_mat(pres, state->pres_prev);
1560 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1562 if ( (membed != nullptr) && (!bLastStep) )
1564 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1567 cycles = wallcycle_stop(wcycle, ewcSTEP);
1568 if (DOMAINDECOMP(cr) && wcycle)
1570 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1573 /* increase the MD step number */
1574 step++;
1575 step_rel++;
1577 resetHandler->resetCounters(
1578 step, step_rel, mdlog, fplog, cr, fr->nbv.get(),
1579 nrnb, fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1581 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1582 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1585 /* End of main MD loop */
1587 /* Closing TNG files can include compressing data. Therefore it is good to do that
1588 * before stopping the time measurements. */
1589 mdoutf_tng_close(outf);
1591 /* Stop measuring walltime */
1592 walltime_accounting_end_time(walltime_accounting);
1594 if (!thisRankHasDuty(cr, DUTY_PME))
1596 /* Tell the PME only node to finish */
1597 gmx_pme_send_finish(cr);
1600 if (MASTER(cr))
1602 if (ir->nstcalcenergy > 0)
1604 energyOutput.printAnnealingTemperatures(fplog, groups, &(ir->opts));
1605 energyOutput.printAverages(fplog, groups);
1608 done_mdoutf(outf);
1610 if (bPMETune)
1612 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1615 done_shellfc(fplog, shellfc, step_rel);
1617 if (useReplicaExchange && MASTER(cr))
1619 print_replica_exchange_statistics(fplog, repl_ex);
1622 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1624 global_stat_destroy(gstat);