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-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/simulation_workload.h"
120 #include "gromacs/mdtypes/state.h"
121 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
122 #include "gromacs/modularsimulator/energyelement.h"
123 #include "gromacs/nbnxm/gpu_data_mgmt.h"
124 #include "gromacs/nbnxm/nbnxm.h"
125 #include "gromacs/pbcutil/mshift.h"
126 #include "gromacs/pbcutil/pbc.h"
127 #include "gromacs/pulling/output.h"
128 #include "gromacs/pulling/pull.h"
129 #include "gromacs/swap/swapcoords.h"
130 #include "gromacs/timing/wallcycle.h"
131 #include "gromacs/timing/walltime_accounting.h"
132 #include "gromacs/topology/atoms.h"
133 #include "gromacs/topology/idef.h"
134 #include "gromacs/topology/mtop_util.h"
135 #include "gromacs/topology/topology.h"
136 #include "gromacs/trajectory/trajectoryframe.h"
137 #include "gromacs/utility/basedefinitions.h"
138 #include "gromacs/utility/cstringutil.h"
139 #include "gromacs/utility/fatalerror.h"
140 #include "gromacs/utility/logger.h"
141 #include "gromacs/utility/real.h"
142 #include "gromacs/utility/smalloc.h"
144 #include "legacysimulator.h"
145 #include "replicaexchange.h"
149 # include "corewrap.h"
152 using gmx::SimulationSignaller;
154 void gmx::LegacySimulator::do_md()
156 // TODO Historically, the EM and MD "integrators" used different
157 // names for the t_inputrec *parameter, but these must have the
158 // same name, now that it's a member of a struct. We use this ir
159 // alias to avoid a large ripple of nearly useless changes.
160 // t_inputrec is being replaced by IMdpOptionsProvider, so this
161 // will go away eventually.
162 t_inputrec* ir = inputrec;
163 int64_t step, step_rel;
164 double t, t0 = ir->init_t, lam0[efptNR];
165 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
166 gmx_bool bNS, bNStList, bStopCM, 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 } }, tmp_vir = { { 0 } },
175 matrix pressureCouplingMu, M;
176 gmx_repl_ex_t repl_ex = nullptr;
178 PaddedHostVector<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 std::vector<RVec> cbuf;
191 real saved_conserved_quantity = 0;
194 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
196 /* PME load balancing data for GPU kernels */
197 gmx_bool bPMETune = FALSE;
198 gmx_bool bPMETunePrinting = FALSE;
200 bool bInteractiveMDstep = false;
202 /* Domain decomposition could incorrectly miss a bonded
203 interaction, but checking for that requires a global
204 communication stage, which does not otherwise happen in DD
205 code. So we do that alongside the first global energy reduction
206 after a new DD is made. These variables handle whether the
207 check happens, and the result it returns. */
208 bool shouldCheckNumberOfBondedInteractions = false;
209 int totalNumberOfBondedInteractions = -1;
211 SimulationSignals signals;
212 // Most global communnication stages don't propagate mdrun
213 // signals, and will use this object to achieve that.
214 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
216 if (!mdrunOptions.writeConfout)
218 // This is on by default, and the main known use case for
219 // turning it off is for convenience in benchmarking, which is
220 // something that should not show up in the general user
225 "The -noconfout functionality is deprecated, and may be removed in a "
229 /* md-vv uses averaged full step velocities for T-control
230 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
231 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
232 bTrotter = (EI_VV(ir->eI)
233 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
235 const bool bRerunMD = false;
237 int nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
238 bGStatEveryStep = (nstglobalcomm == 1);
240 SimulationGroups* groups = &top_global->groups;
242 std::unique_ptr<EssentialDynamics> ed = nullptr;
243 if (opt2bSet("-ei", nfile, fnm))
245 /* Initialize essential dynamics sampling */
246 ed = init_edsam(mdlog, opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm), top_global,
247 ir, cr, constr, state_global, observablesHistory, oenv, startingBehavior);
249 else if (observablesHistory->edsamHistory)
252 "The checkpoint is from a run with essential dynamics sampling, "
253 "but the current run did not specify the -ei option. "
254 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
257 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
258 Update upd(ir, deform);
259 const bool doSimulatedAnnealing = initSimulatedAnnealing(ir, &upd);
260 if (startingBehavior != StartingBehavior::RestartWithAppending)
262 pleaseCiteCouplingAlgorithms(fplog, *ir);
264 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
265 mdModulesNotifier, ir, top_global, oenv, wcycle, startingBehavior);
266 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
267 mdoutf_get_fp_dhdl(outf), false, startingBehavior, mdModulesNotifier);
269 gstat = global_stat_init(ir);
271 /* Check for polarizable models and flexible constraints */
272 shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
273 ir->nstcalcenergy, DOMAINDECOMP(cr));
276 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
277 if ((io > 2000) && MASTER(cr))
279 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
283 // Local state only becomes valid now.
284 std::unique_ptr<t_state> stateInstance;
288 auto mdatoms = mdAtoms->mdatoms();
290 std::unique_ptr<UpdateConstrainCuda> integrator;
292 if (DOMAINDECOMP(cr))
294 dd_init_local_top(*top_global, &top);
296 stateInstance = std::make_unique<t_state>();
297 state = stateInstance.get();
298 dd_init_local_state(cr->dd, state_global, state);
300 /* Distribute the charge groups over the nodes from the master node */
301 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
302 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
303 nrnb, nullptr, FALSE);
304 shouldCheckNumberOfBondedInteractions = true;
305 upd.setNumAtoms(state->natoms);
309 state_change_natoms(state_global, state_global->natoms);
310 f.resizeWithPadding(state_global->natoms);
311 /* Copy the pointer to the global state */
312 state = state_global;
314 /* Generate and initialize new topology */
315 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &graph, mdAtoms, constr, vsite, shellfc);
317 upd.setNumAtoms(state->natoms);
320 const auto& simulationWork = runScheduleWork->simulationWork;
321 const bool useGpuForPme = simulationWork.useGpuPme;
322 const bool useGpuForNonbonded = simulationWork.useGpuNonbonded;
323 const bool useGpuForBufferOps = simulationWork.useGpuBufferOps;
324 const bool useGpuForUpdate = simulationWork.useGpuUpdate;
326 StatePropagatorDataGpu* stateGpu = fr->stateGpu;
330 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
331 || constr->numConstraintsTotal() == 0,
332 "Constraints in domain decomposition are only supported with update "
333 "groups if using GPU update.\n");
334 GMX_RELEASE_ASSERT(ir->eConstrAlg != econtSHAKE || constr == nullptr
335 || constr->numConstraintsTotal() == 0,
336 "SHAKE is not supported with GPU update.");
337 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
338 "Either PME or short-ranged non-bonded interaction tasks must run on "
339 "the GPU to use GPU update.\n");
340 GMX_RELEASE_ASSERT(ir->eI == eiMD,
341 "Only the md integrator is supported with the GPU update.\n");
343 ir->etc != etcNOSEHOOVER,
344 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
345 GMX_RELEASE_ASSERT(ir->epc == epcNO || ir->epc == epcPARRINELLORAHMAN || ir->epc == epcBERENDSEN,
346 "Only Parrinello-Rahman and Berendsen pressure coupling are supported "
347 "with the GPU update.\n");
348 GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
349 "Virtual sites are not supported with the GPU update.\n");
350 GMX_RELEASE_ASSERT(ed == nullptr,
351 "Essential dynamics is not supported with the GPU update.\n");
352 GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(ir->pull),
353 "Constraints pulling is not supported with the GPU update.\n");
354 GMX_RELEASE_ASSERT(fcd->orires.nr == 0,
355 "Orientation restraints are not supported with the GPU update.\n");
356 GMX_RELEASE_ASSERT(ir->efep == efepNO,
357 "Free energy perturbations are not supported with the GPU update.");
358 GMX_RELEASE_ASSERT(graph == nullptr, "The graph is not supported with GPU update.");
360 if (constr != nullptr && constr->numConstraintsTotal() > 0)
364 .appendText("Updating coordinates and applying constraints on the GPU.");
368 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
370 integrator = std::make_unique<UpdateConstrainCuda>(
371 *ir, *top_global, stateGpu->getUpdateStream(), stateGpu->xUpdatedOnDevice());
374 set_pbc(&pbc, epbcXYZ, state->box);
375 integrator->setPbc(&pbc);
378 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
380 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
382 if ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
384 changePinningPolicy(&f, PinningPolicy::PinnedIfSupported);
388 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
391 // NOTE: The global state is no longer used at this point.
392 // But state_global is still used as temporary storage space for writing
393 // the global state to file and potentially for replica exchange.
394 // (Global topology should persist.)
396 update_mdatoms(mdatoms, state->lambda[efptMASS]);
400 /* Check nstexpanded here, because the grompp check was broken */
401 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
404 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
406 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
411 EnergyElement::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
414 preparePrevStepPullCom(ir, pull_work, mdatoms, state, state_global, cr,
415 startingBehavior != StartingBehavior::NewSimulation);
417 // TODO: Remove this by converting AWH into a ForceProvider
418 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms,
419 startingBehavior != StartingBehavior::NewSimulation,
420 shellfc != nullptr, opt2fn("-awh", nfile, fnm), pull_work);
422 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
423 if (useReplicaExchange && MASTER(cr))
425 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir, replExParams);
427 /* PME tuning is only supported in the Verlet scheme, with PME for
428 * Coulomb. It is not supported with only LJ PME. */
429 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
430 && ir->cutoff_scheme != ecutsGROUP);
432 pme_load_balancing_t* pme_loadbal = nullptr;
435 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata,
439 if (!ir->bContinuation)
441 if (state->flags & (1U << estV))
443 auto v = makeArrayRef(state->v);
444 /* Set the velocities of vsites, shells and frozen atoms to zero */
445 for (i = 0; i < mdatoms->homenr; i++)
447 if (mdatoms->ptype[i] == eptVSite || mdatoms->ptype[i] == eptShell)
451 else if (mdatoms->cFREEZE)
453 for (m = 0; m < DIM; m++)
455 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
466 /* Constrain the initial coordinates and velocities */
467 do_constrain_first(fplog, constr, ir, mdatoms, state->natoms, state->x.arrayRefWithPadding(),
468 state->v.arrayRefWithPadding(), state->box, state->lambda[efptBONDED]);
472 /* Construct the virtual sites for the initial configuration */
473 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, nullptr, top.idef.iparams,
474 top.idef.il, fr->ePBC, fr->bMolPBC, cr, state->box);
478 if (ir->efep != efepNO)
480 /* Set free energy calculation frequency as the greatest common
481 * denominator of nstdhdl and repl_ex_nst. */
482 nstfep = ir->fepvals->nstdhdl;
485 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
487 if (useReplicaExchange)
489 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
493 /* Be REALLY careful about what flags you set here. You CANNOT assume
494 * this is the first step, since we might be restarting from a checkpoint,
495 * and in that case we should not do any modifications to the state.
497 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
499 // When restarting from a checkpoint, it can be appropriate to
500 // initialize ekind from quantities in the checkpoint. Otherwise,
501 // compute_globals must initialize ekind before the simulation
502 // starts/restarts. However, only the master rank knows what was
503 // found in the checkpoint file, so we have to communicate in
504 // order to coordinate the restart.
506 // TODO Consider removing this communication if/when checkpoint
507 // reading directly follows .tpr reading, because all ranks can
508 // agree on hasReadEkinState at that time.
509 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
512 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr);
514 if (hasReadEkinState)
516 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
519 unsigned int cglo_flags =
520 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
521 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (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, state->x.rvec_array(), state->v.rvec_array(),
542 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, nullptr, enerd,
543 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
544 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
546 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
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, top_global, &top,
563 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, state->x.rvec_array(), state->v.rvec_array(),
574 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, nullptr, enerd,
575 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
576 state->box, nullptr, &bSumEkinhOld, cglo_flags & ~CGLO_PRESSURE);
579 /* Calculate the initial half step temperature, and save the ekinh_old */
580 if (startingBehavior == StartingBehavior::NewSimulation)
582 for (i = 0; (i < ir->opts.ngtc); i++)
584 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
588 /* need to make an initiation call to get the Trotter variables set, as well as other constants
589 for non-trotter temperature control */
590 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
594 if (!ir->bContinuation)
596 if (constr && ir->eConstrAlg == econtLINCS)
598 fprintf(fplog, "RMS relative constraint deviation after constraining: %.2e\n",
601 if (EI_STATE_VELOCITY(ir->eI))
603 real temp = enerd->term[F_TEMP];
606 /* Result of Ekin averaged over velocities of -half
607 * and +half step, while we only have -half step here.
611 fprintf(fplog, "Initial temperature: %g K\n", temp);
616 fprintf(stderr, "starting mdrun '%s'\n", *(top_global->name));
619 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
623 sprintf(tbuf, "%s", "infinite");
625 if (ir->init_step > 0)
627 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
628 gmx_step_str(ir->init_step + ir->nsteps, sbuf), tbuf,
629 gmx_step_str(ir->init_step, sbuf2), ir->init_step * ir->delta_t);
633 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
635 fprintf(fplog, "\n");
638 walltime_accounting_start_time(walltime_accounting);
639 wallcycle_start(wcycle, ewcRUN);
640 print_start(fplog, cr, walltime_accounting, "mdrun");
643 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
644 int chkpt_ret = fcCheckPointParallel(cr->nodeid, NULL, 0);
647 gmx_fatal(3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0);
651 /***********************************************************
655 ************************************************************/
658 /* Skip the first Nose-Hoover integration when we get the state from tpx */
659 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
660 bSumEkinhOld = FALSE;
662 bNeedRepartition = FALSE;
664 bool simulationsShareState = false;
665 int nstSignalComm = nstglobalcomm;
667 // TODO This implementation of ensemble orientation restraints is nasty because
668 // a user can't just do multi-sim with single-sim orientation restraints.
669 bool usingEnsembleRestraints =
670 (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
671 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
673 // Replica exchange, ensemble restraints and AWH need all
674 // simulations to remain synchronized, so they need
675 // checkpoints and stop conditions to act on the same step, so
676 // the propagation of such signals must take place between
677 // simulations, not just within simulations.
678 // TODO: Make algorithm initializers set these flags.
679 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
681 if (simulationsShareState)
683 // Inter-simulation signal communication does not need to happen
684 // often, so we use a minimum of 200 steps to reduce overhead.
685 const int c_minimumInterSimulationSignallingInterval = 200;
686 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
691 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
692 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
693 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
694 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
696 auto checkpointHandler = std::make_unique<CheckpointHandler>(
697 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]), simulationsShareState,
698 ir->nstlist == 0, MASTER(cr), mdrunOptions.writeConfout,
699 mdrunOptions.checkpointOptions.period);
701 const bool resetCountersIsLocal = true;
702 auto resetHandler = std::make_unique<ResetHandler>(
703 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
704 !resetCountersIsLocal, ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
705 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
707 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
709 step = ir->init_step;
712 // TODO extract this to new multi-simulation module
713 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
715 if (!multisim_int_all_are_equal(ms, ir->nsteps))
717 GMX_LOG(mdlog.warning)
719 "Note: The number of steps is not consistent across multi "
721 "but we are proceeding anyway!");
723 if (!multisim_int_all_are_equal(ms, ir->init_step))
725 GMX_LOG(mdlog.warning)
727 "Note: The initial step is not consistent across multi simulations,\n"
728 "but we are proceeding anyway!");
732 /* and stop now if we should */
733 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
737 /* Determine if this is a neighbor search step */
738 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
740 if (bPMETune && bNStList)
742 // This has to be here because PME load balancing is called so early.
743 // TODO: Move to after all booleans are defined.
744 if (useGpuForUpdate && !bFirstStep)
746 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
747 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
749 /* PME grid + cut-off optimization with GPUs or PME nodes */
750 pme_loadbal_do(pme_loadbal, cr, (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
751 fplog, mdlog, *ir, fr, state->box, state->x, wcycle, step, step_rel,
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->lambda, state->fep_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) && (ir->bExpanded)
769 && (step > 0) && (startingBehavior == StartingBehavior::NewSimulation));
772 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
773 && do_per_step(step, replExParams.exchangeInterval));
775 if (doSimulatedAnnealing)
777 update_annealing_target_temp(ir, t, &upd);
780 /* Stop Center of Mass motion */
781 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
783 /* Determine whether or not to do Neighbour Searching */
784 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
786 /* Note that the stopHandler will cause termination at nstglobalcomm
787 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
788 * nstpcouple steps, we have computed the half-step kinetic energy
789 * of the previous step and can always output energies at the last step.
791 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
793 /* do_log triggers energy and virial calculation. Because this leads
794 * to different code paths, forces can be different. Thus for exact
795 * continuation we should avoid extra log output.
796 * Note that the || bLastStep can result in non-exact continuation
797 * beyond the last step. But we don't consider that to be an issue.
799 do_log = (do_per_step(step, ir->nstlog)
800 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
801 do_verbose = mdrunOptions.verbose
802 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
804 if (useGpuForUpdate && !bFirstStep && bNS)
806 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
807 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
808 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
809 // Copy coordinate from the GPU when needed at the search step.
810 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
811 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
812 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
813 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
816 if (bNS && !(bFirstStep && ir->bContinuation))
818 bMasterState = FALSE;
819 /* Correct the new box if it is too skewed */
820 if (inputrecDynamicBox(ir))
822 if (correct_box(fplog, step, state->box, graph))
825 // If update is offloaded, it should be informed about the box size change
829 set_pbc(&pbc, epbcXYZ, state->box);
830 integrator->setPbc(&pbc);
834 if (DOMAINDECOMP(cr) && bMasterState)
836 dd_collect_state(cr->dd, state, state_global);
839 if (DOMAINDECOMP(cr))
841 /* Repartition the domain decomposition */
842 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
843 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
844 fr, vsite, constr, nrnb, wcycle, do_verbose && !bPMETunePrinting);
845 shouldCheckNumberOfBondedInteractions = true;
846 upd.setNumAtoms(state->natoms);
850 if (MASTER(cr) && do_log)
852 energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
855 if (ir->efep != efepNO)
857 update_mdatoms(mdatoms, state->lambda[efptMASS]);
863 /* We need the kinetic energy at minus the half step for determining
864 * the full step kinetic energy and possibly for T-coupling.*/
865 /* This may not be quite working correctly yet . . . . */
866 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
867 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
868 nullptr, nullptr, nullptr, nullptr, mu_tot, constr, &nullSignaller,
869 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
870 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
871 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
872 &top, state->x.rvec_array(), state->box,
873 &shouldCheckNumberOfBondedInteractions);
875 clear_mat(force_vir);
877 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
879 /* Determine the energy and pressure:
880 * at nstcalcenergy steps and at energy output steps (set below).
882 if (EI_VV(ir->eI) && (!bInitStep))
884 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
885 bCalcVir = bCalcEnerStep
887 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
891 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
892 bCalcVir = bCalcEnerStep || (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
894 bCalcEner = bCalcEnerStep;
896 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
898 if (do_ene || do_log || bDoReplEx)
904 /* Do we need global communication ? */
905 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
906 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
908 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
909 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
910 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
914 /* Now is the time to relax the shells */
915 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
916 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
917 state->natoms, state->x.arrayRefWithPadding(),
918 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
919 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, graph,
920 shellfc, fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
924 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
925 is updated (or the AWH update will be performed twice for one step when continuing).
926 It would be best to call this update function from do_md_trajectory_writing but that
927 would occur after do_force. One would have to divide the update_awh function into one
928 function applying the AWH force and one doing the AWH bias update. The update AWH
929 bias function could then be called after do_md_trajectory_writing (then containing
930 update_awh_history). The checkpointing will in the future probably moved to the start
931 of the md loop which will 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)
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, pull_work, step,
943 nrnb, wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
944 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, graph,
945 fr, runScheduleWork, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
946 (bNS ? GMX_FORCE_NS : 0) | force_flags, ddBalanceRegionHandler);
949 // VV integrators do not need the following velocity half step
950 // if it is the first step after starting from a checkpoint.
951 // That is, the half step is needed on all other steps, and
952 // also the first step when starting from a .tpr file.
953 if (EI_VV(ir->eI) && (!bFirstStep || 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,
969 state->natoms); /* should make this better for parallelizing? */
973 /* this is for NHC in the Ekin(t+dt/2) version of vv */
974 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
975 trotter_seq, ettTSEQ1);
978 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
979 etrtVELOCITY1, cr, constr);
981 wallcycle_stop(wcycle, ewcUPDATE);
982 constrain_velocities(step, nullptr, state, shake_vir, constr, bCalcVir, do_log, do_ene);
983 wallcycle_start(wcycle, ewcUPDATE);
984 /* if VV, compute the pressure and constraints */
985 /* For VV2, we strictly only need this if using pressure
986 * control, but we really would like to have accurate pressures
988 * Think about ways around this in the future?
989 * For now, keep this choice in comments.
991 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
992 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
994 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
995 if (bCalcEner && ir->eI == eiVVAK)
999 /* for vv, the first half of the integration actually corresponds to the previous step.
1000 So we need information from the last step in the first half of the integration */
1001 if (bGStat || do_per_step(step - 1, nstglobalcomm))
1003 wallcycle_stop(wcycle, ewcUPDATE);
1004 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1005 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1006 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
1007 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1008 (bGStat ? CGLO_GSTAT : 0) | (bCalcEner ? CGLO_ENERGY : 0)
1009 | (bTemp ? CGLO_TEMPERATURE : 0) | (bPres ? CGLO_PRESSURE : 0)
1010 | (bPres ? CGLO_CONSTRAINT : 0) | (bStopCM ? CGLO_STOPCM : 0)
1011 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1014 /* explanation of above:
1015 a) We compute Ekin at the full time step
1016 if 1) we are using the AveVel Ekin, and it's not the
1017 initial step, or 2) if we are using AveEkin, but need the full
1018 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1019 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1020 EkinAveVel because it's needed for the pressure */
1021 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1022 top_global, &top, state->x.rvec_array(), state->box,
1023 &shouldCheckNumberOfBondedInteractions);
1026 process_and_stopcm_grp(fplog, &vcm, *mdatoms, state->x.rvec_array(),
1027 state->v.rvec_array());
1028 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1030 wallcycle_start(wcycle, ewcUPDATE);
1032 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1037 m_add(force_vir, shake_vir,
1038 total_vir); /* we need the un-dispersion corrected total vir here */
1039 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1040 trotter_seq, ettTSEQ2);
1042 /* TODO This is only needed when we're about to write
1043 * a checkpoint, because we use it after the restart
1044 * (in a kludge?). But what should we be doing if
1045 * the startingBehavior is NewSimulation or bInitStep are true? */
1046 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1048 copy_mat(shake_vir, state->svir_prev);
1049 copy_mat(force_vir, state->fvir_prev);
1051 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1053 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1054 enerd->term[F_TEMP] =
1055 sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1056 enerd->term[F_EKIN] = trace(ekind->ekin);
1059 else if (bExchanged)
1061 wallcycle_stop(wcycle, ewcUPDATE);
1062 /* We need the kinetic energy at minus the half step for determining
1063 * the full step kinetic energy and possibly for T-coupling.*/
1064 /* This may not be quite working correctly yet . . . . */
1065 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(),
1066 state->v.rvec_array(), state->box, state->lambda[efptVDW],
1067 mdatoms, nrnb, &vcm, wcycle, enerd, nullptr, nullptr, nullptr,
1068 nullptr, mu_tot, constr, &nullSignaller, state->box, nullptr,
1069 &bSumEkinhOld, CGLO_GSTAT | CGLO_TEMPERATURE);
1070 wallcycle_start(wcycle, ewcUPDATE);
1073 /* if it's the initial step, we performed this first step just to get the constraint virial */
1074 if (ir->eI == eiVV && bInitStep)
1076 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1079 wallcycle_stop(wcycle, ewcUPDATE);
1082 /* compute the conserved quantity */
1085 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1088 last_ekin = enerd->term[F_EKIN];
1090 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1092 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1094 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1095 if (ir->efep != efepNO)
1097 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1101 /* ######## END FIRST UPDATE STEP ############## */
1102 /* ######## If doing VV, we now have v(dt) ###### */
1105 /* perform extended ensemble sampling in lambda - we don't
1106 actually move to the new state before outputting
1107 statistics, but if performing simulated tempering, we
1108 do update the velocities and the tau_t. */
1110 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state,
1111 state->dfhist, step, state->v.rvec_array(), mdatoms);
1112 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1115 copy_df_history(state_global->dfhist, state->dfhist);
1119 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1120 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1121 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1122 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1123 || checkpointHandler->isCheckpointingStep()))
1125 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1126 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1128 // Copy velocities if needed for the output/checkpointing.
1129 // NOTE: Copy on the search steps is done at the beginning of the step.
1130 if (useGpuForUpdate && !bNS
1131 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1133 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1134 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1136 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1137 // and update is offloaded hence forces are kept on the GPU for update and have not been
1138 // already transferred in do_force().
1139 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1140 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1141 // prior to GPU update.
1142 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1143 // copy call in do_force(...).
1144 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1145 // on host after the D2H copy in do_force(...).
1146 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1147 && do_per_step(step, ir->nstfout))
1149 stateGpu->copyForcesFromGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1150 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1152 /* Now we have the energies and forces corresponding to the
1153 * coordinates at time t. We must output all of this before
1156 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state, state_global,
1157 observablesHistory, top_global, fr, outf, energyOutput, ekind, f,
1158 checkpointHandler->isCheckpointingStep(), bRerunMD, bLastStep,
1159 mdrunOptions.writeConfout, bSumEkinhOld);
1160 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1161 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1163 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1164 if (startingBehavior != StartingBehavior::NewSimulation
1165 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1167 copy_mat(state->svir_prev, shake_vir);
1168 copy_mat(state->fvir_prev, force_vir);
1171 stopHandler->setSignal();
1172 resetHandler->setSignal(walltime_accounting);
1174 if (bGStat || !PAR(cr))
1176 /* In parallel we only have to check for checkpointing in steps
1177 * where we do global communication,
1178 * otherwise the other nodes don't know.
1180 checkpointHandler->setSignal(walltime_accounting);
1183 /* ######### START SECOND UPDATE STEP ################# */
1185 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1186 controlled in preprocessing */
1188 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1190 gmx_bool bIfRandomize;
1191 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1192 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1193 if (constr && bIfRandomize)
1195 constrain_velocities(step, nullptr, state, tmp_vir, constr, bCalcVir, do_log, do_ene);
1198 /* Box is changed in update() when we do pressure coupling,
1199 * but we should still use the old box for energy corrections and when
1200 * writing it to the energy file, so it matches the trajectory files for
1201 * the same timestep above. Make a copy in a separate array.
1203 copy_mat(state->box, lastbox);
1207 wallcycle_start(wcycle, ewcUPDATE);
1208 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1211 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1212 /* We can only do Berendsen coupling after we have summed
1213 * the kinetic energy or virial. Since the happens
1214 * in global_state after update, we should only do it at
1215 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1220 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1221 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1226 /* velocity half-step update */
1227 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1228 etrtVELOCITY2, cr, constr);
1231 /* Above, initialize just copies ekinh into ekin,
1232 * it doesn't copy position (for VV),
1233 * and entire integrator for MD.
1236 if (ir->eI == eiVVAK)
1238 cbuf.resize(state->x.size());
1239 std::copy(state->x.begin(), state->x.end(), cbuf.begin());
1242 /* With leap-frog type integrators we compute the kinetic energy
1243 * at a whole time step as the average of the half-time step kinetic
1244 * energies of two subsequent steps. Therefore we need to compute the
1245 * half step kinetic energy also if we need energies at the next step.
1247 const bool needHalfStepKineticEnergy =
1248 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1250 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1251 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1252 const bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN
1253 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1255 if (useGpuForUpdate)
1259 integrator->set(stateGpu->getCoordinates(), stateGpu->getVelocities(),
1260 stateGpu->getForces(), top.idef, *mdatoms, ekind->ngtc);
1262 // Copy data to the GPU after buffers might have being reinitialized
1263 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1264 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1267 // If the buffer ops were not offloaded this step, the forces are on the host and have to be copied
1268 if (!runScheduleWork->stepWork.useGpuFBufferOps)
1270 stateGpu->copyForcesToGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1273 const bool doTemperatureScaling =
1274 (ir->etc != etcNO && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1276 // This applies Leap-Frog, LINCS and SETTLE in succession
1277 integrator->integrate(stateGpu->getForcesReadyOnDeviceEvent(
1278 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1279 ir->delta_t, true, bCalcVir, shake_vir, doTemperatureScaling,
1280 ekind->tcstat, doParrinelloRahman, ir->nstpcouple * ir->delta_t, M);
1282 // Copy velocities D2H after update if:
1283 // - Globals are computed this step (includes the energy output steps).
1284 // - Temperature is needed for the next step.
1285 if (bGStat || needHalfStepKineticEnergy)
1287 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1288 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1293 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1294 etrtPOSITION, cr, constr);
1296 wallcycle_stop(wcycle, ewcUPDATE);
1298 constrain_coordinates(step, &dvdl_constr, state, shake_vir, &upd, constr, bCalcVir,
1301 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state, cr, nrnb, wcycle, &upd,
1302 constr, do_log, do_ene);
1303 finish_update(ir, mdatoms, state, graph, nrnb, wcycle, &upd, constr);
1306 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1308 updatePrevStepPullCom(pull_work, state);
1311 if (ir->eI == eiVVAK)
1313 /* erase F_EKIN and F_TEMP here? */
1314 /* just compute the kinetic energy at the half step to perform a trotter step */
1315 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1316 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1317 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller, lastbox,
1318 nullptr, &bSumEkinhOld, (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE);
1319 wallcycle_start(wcycle, ewcUPDATE);
1320 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1321 /* now we know the scaling, we can compute the positions again */
1322 std::copy(cbuf.begin(), cbuf.end(), state->x.begin());
1324 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1325 etrtPOSITION, cr, constr);
1326 wallcycle_stop(wcycle, ewcUPDATE);
1328 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1329 /* are the small terms in the shake_vir here due
1330 * to numerical errors, or are they important
1331 * physically? I'm thinking they are just errors, but not completely sure.
1332 * For now, will call without actually constraining, constr=NULL*/
1333 finish_update(ir, mdatoms, state, graph, nrnb, wcycle, &upd, nullptr);
1337 /* this factor or 2 correction is necessary
1338 because half of the constraint force is removed
1339 in the vv step, so we have to double it. See
1340 the Redmine issue #1255. It is not yet clear
1341 if the factor of 2 is exact, or just a very
1342 good approximation, and this will be
1343 investigated. The next step is to see if this
1344 can be done adding a dhdl contribution from the
1345 rattle step, but this is somewhat more
1346 complicated with the current code. Will be
1347 investigated, hopefully for 4.6.3. However,
1348 this current solution is much better than
1349 having it completely wrong.
1351 enerd->term[F_DVDL_CONSTR] += 2 * dvdl_constr;
1355 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1358 if (vsite != nullptr)
1360 wallcycle_start(wcycle, ewcVSITECONSTR);
1361 if (graph != nullptr)
1363 shift_self(graph, state->box, state->x.rvec_array());
1365 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, state->v.rvec_array(),
1366 top.idef.iparams, top.idef.il, fr->ePBC, fr->bMolPBC, cr, state->box);
1368 if (graph != nullptr)
1370 unshift_self(graph, state->box, state->x.rvec_array());
1372 wallcycle_stop(wcycle, ewcVSITECONSTR);
1375 /* ############## IF NOT VV, Calculate globals HERE ############ */
1376 /* With Leap-Frog we can skip compute_globals at
1377 * non-communication steps, but we need to calculate
1378 * the kinetic energy one step before communication.
1381 // Organize to do inter-simulation signalling on steps if
1382 // and when algorithms require it.
1383 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1385 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1387 // Copy coordinates when needed to stop the CM motion.
1388 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1390 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1391 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1393 // Since we're already communicating at this step, we
1394 // can propagate intra-simulation signals. Note that
1395 // check_nstglobalcomm has the responsibility for
1396 // choosing the value of nstglobalcomm that is one way
1397 // bGStat becomes true, so we can't get into a
1398 // situation where e.g. checkpointing can't be
1400 bool doIntraSimSignal = true;
1401 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1404 gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1405 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1406 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &signaller, lastbox,
1407 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1408 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1409 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1410 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1411 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1412 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1414 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1415 top_global, &top, state->x.rvec_array(), state->box,
1416 &shouldCheckNumberOfBondedInteractions);
1417 if (!EI_VV(ir->eI) && bStopCM)
1419 process_and_stopcm_grp(fplog, &vcm, *mdatoms, state->x.rvec_array(),
1420 state->v.rvec_array());
1421 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1423 // TODO: The special case of removing CM motion should be dealt more gracefully
1424 if (useGpuForUpdate)
1426 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1427 // Here we block until the H2D copy completes because event sync with the
1428 // force kernels that use the coordinates on the next steps is not implemented
1429 // (not because of a race on state->x being modified on the CPU while H2D is in progress).
1430 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1436 /* ############# END CALC EKIN AND PRESSURE ################# */
1438 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1439 the virial that should probably be addressed eventually. state->veta has better properies,
1440 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1441 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1443 if (ir->efep != efepNO && !EI_VV(ir->eI))
1445 /* Sum up the foreign energy and dhdl terms for md and sd.
1446 Currently done every step so that dhdl is correct in the .edr */
1447 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1450 update_pcouple_after_coordinates(fplog, step, ir, mdatoms, pres, force_vir, shake_vir,
1451 pressureCouplingMu, state, nrnb, &upd, !useGpuForUpdate);
1453 const bool doBerendsenPressureCoupling =
1454 (inputrec->epc == epcBERENDSEN && do_per_step(step, inputrec->nstpcouple));
1455 if (useGpuForUpdate && (doBerendsenPressureCoupling || doParrinelloRahman))
1457 integrator->scaleCoordinates(pressureCouplingMu);
1459 set_pbc(&pbc, epbcXYZ, state->box);
1460 integrator->setPbc(&pbc);
1463 /* ################# END UPDATE STEP 2 ################# */
1464 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1466 /* The coordinates (x) were unshifted in update */
1469 /* We will not sum ekinh_old,
1470 * so signal that we still have to do it.
1472 bSumEkinhOld = TRUE;
1477 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1479 /* use the directly determined last velocity, not actually the averaged half steps */
1480 if (bTrotter && ir->eI == eiVV)
1482 enerd->term[F_EKIN] = last_ekin;
1484 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1486 if (integratorHasConservedEnergyQuantity(ir))
1490 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1494 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1497 /* ######### END PREPARING EDR OUTPUT ########### */
1503 if (fplog && do_log && bDoExpanded)
1505 /* only needed if doing expanded ensemble */
1506 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals,
1507 ir->bSimTemp ? ir->simtempvals : nullptr,
1508 state_global->dfhist, state->fep_state, ir->nstlog, step);
1512 energyOutput.addDataAtEnergyStep(bDoDHDL, bCalcEnerStep, t, mdatoms->tmass, enerd, state,
1513 ir->fepvals, ir->expandedvals, lastbox, shake_vir,
1514 force_vir, total_vir, pres, ekind, mu_tot, constr);
1518 energyOutput.recordNonEnergyStep();
1521 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1522 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1524 if (doSimulatedAnnealing)
1526 energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
1528 if (do_log || do_ene || do_dr || do_or)
1530 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1531 do_log ? fplog : nullptr, step, t, fcd, awh.get());
1536 pull_print_output(pull_work, step, t);
1539 if (do_per_step(step, ir->nstlog))
1541 if (fflush(fplog) != 0)
1543 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1549 /* Have to do this part _after_ outputting the logfile and the edr file */
1550 /* Gets written into the state at the beginning of next loop*/
1551 state->fep_state = lamnew;
1553 /* Print the remaining wall clock time for the run */
1554 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1558 fprintf(stderr, "\n");
1560 print_time(stderr, walltime_accounting, step, ir, cr);
1563 /* Ion/water position swapping.
1564 * Not done in last step since trajectory writing happens before this call
1565 * in the MD loop and exchanges would be lost anyway. */
1566 bNeedRepartition = FALSE;
1567 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep && do_per_step(step, ir->swap->nstswap))
1570 do_swapcoords(cr, step, t, ir, swap, wcycle, as_rvec_array(state->x.data()),
1571 state->box, MASTER(cr) && mdrunOptions.verbose, bRerunMD);
1573 if (bNeedRepartition && DOMAINDECOMP(cr))
1575 dd_collect_state(cr->dd, state, state_global);
1579 /* Replica exchange */
1583 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1586 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1588 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1, state_global, *top_global, ir,
1589 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
1590 nrnb, wcycle, FALSE);
1591 shouldCheckNumberOfBondedInteractions = true;
1592 upd.setNumAtoms(state->natoms);
1598 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1599 /* With all integrators, except VV, we need to retain the pressure
1600 * at the current step for coupling at the next step.
1602 if ((state->flags & (1U << estPRES_PREV))
1603 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1605 /* Store the pressure in t_state for pressure coupling
1606 * at the next MD step.
1608 copy_mat(pres, state->pres_prev);
1611 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1613 if ((membed != nullptr) && (!bLastStep))
1615 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1618 cycles = wallcycle_stop(wcycle, ewcSTEP);
1619 if (DOMAINDECOMP(cr) && wcycle)
1621 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1624 /* increase the MD step number */
1628 resetHandler->resetCounters(step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb,
1629 fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1631 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1632 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1634 /* End of main MD loop */
1636 /* Closing TNG files can include compressing data. Therefore it is good to do that
1637 * before stopping the time measurements. */
1638 mdoutf_tng_close(outf);
1640 /* Stop measuring walltime */
1641 walltime_accounting_end_time(walltime_accounting);
1643 if (!thisRankHasDuty(cr, DUTY_PME))
1645 /* Tell the PME only node to finish */
1646 gmx_pme_send_finish(cr);
1651 if (ir->nstcalcenergy > 0)
1653 energyOutput.printAnnealingTemperatures(fplog, groups, &(ir->opts));
1654 energyOutput.printAverages(fplog, groups);
1661 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1664 done_shellfc(fplog, shellfc, step_rel);
1666 if (useReplicaExchange && MASTER(cr))
1668 print_replica_exchange_statistics(fplog, repl_ex);
1671 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1673 global_stat_destroy(gstat);