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, 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 /*****************************************************************************************/
321 // TODO: The following block of code should be refactored, once:
322 // 1. We have the useGpuForBufferOps variable set and available here and in do_force(...)
323 // 2. The proper GPU syncronization is introduced, so that the H2D and D2H data copies can be performed in the separate
324 // stream owned by the StatePropagatorDataGpu
325 const auto& simulationWork = runScheduleWork->simulationWork;
326 const bool useGpuForPme = simulationWork.useGpuPme;
327 const bool useGpuForNonbonded = simulationWork.useGpuNonbonded;
328 // Temporary solution to make sure that the buffer ops are offloaded when update is offloaded
329 const bool useGpuForBufferOps = simulationWork.useGpuBufferOps;
330 const bool useGpuForUpdate = simulationWork.useGpuUpdate;
333 StatePropagatorDataGpu* stateGpu = fr->stateGpu;
337 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || ddUsesUpdateGroups(*cr->dd) || constr == nullptr
338 || constr->numConstraintsTotal() == 0,
339 "Constraints in domain decomposition are only supported with update "
340 "groups if using GPU update.\n");
341 GMX_RELEASE_ASSERT(ir->eConstrAlg != econtSHAKE || constr == nullptr
342 || constr->numConstraintsTotal() == 0,
343 "SHAKE is not supported with GPU update.");
344 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
345 "Either PME or short-ranged non-bonded interaction tasks must run on "
346 "the GPU to use GPU update.\n");
347 GMX_RELEASE_ASSERT(ir->eI == eiMD,
348 "Only the md integrator is supported with the GPU update.\n");
350 ir->etc != etcNOSEHOOVER,
351 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
352 GMX_RELEASE_ASSERT(ir->epc == epcNO || ir->epc == epcPARRINELLORAHMAN || ir->epc == epcBERENDSEN,
353 "Only Parrinello-Rahman and Berendsen pressure coupling are supported "
354 "with the GPU update.\n");
355 GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
356 "Virtual sites are not supported with the GPU update.\n");
357 GMX_RELEASE_ASSERT(ed == nullptr,
358 "Essential dynamics is not supported with the GPU update.\n");
359 GMX_RELEASE_ASSERT(!ir->bPull || !pull_have_constraint(ir->pull),
360 "Pulling is not supported with the GPU update.\n");
361 GMX_RELEASE_ASSERT(fcd->orires.nr == 0,
362 "Orientation restraints are not supported with the GPU update.\n");
363 GMX_RELEASE_ASSERT(ir->efep == efepNO,
364 "Free energy perturbations are not supported with the GPU update.");
365 GMX_RELEASE_ASSERT(graph == nullptr, "The graph is not supported with GPU update.");
367 if (constr != nullptr && constr->numConstraintsTotal() > 0)
371 .appendText("Updating coordinates and applying constraints on the GPU.");
375 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
377 integrator = std::make_unique<UpdateConstrainCuda>(
378 *ir, *top_global, stateGpu->getUpdateStream(), stateGpu->xUpdatedOnDevice());
381 set_pbc(&pbc, epbcXYZ, state->box);
382 integrator->setPbc(&pbc);
385 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
387 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
389 if ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
391 changePinningPolicy(&f, PinningPolicy::PinnedIfSupported);
395 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
397 /*****************************************************************************************/
399 // NOTE: The global state is no longer used at this point.
400 // But state_global is still used as temporary storage space for writing
401 // the global state to file and potentially for replica exchange.
402 // (Global topology should persist.)
404 update_mdatoms(mdatoms, state->lambda[efptMASS]);
408 /* Check nstexpanded here, because the grompp check was broken */
409 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
412 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
414 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
419 EnergyElement::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
422 preparePrevStepPullCom(ir, pull_work, mdatoms, state, state_global, cr,
423 startingBehavior != StartingBehavior::NewSimulation);
425 // TODO: Remove this by converting AWH into a ForceProvider
426 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms,
427 startingBehavior != StartingBehavior::NewSimulation,
428 shellfc != nullptr, opt2fn("-awh", nfile, fnm), pull_work);
430 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
431 if (useReplicaExchange && MASTER(cr))
433 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir, replExParams);
435 /* PME tuning is only supported in the Verlet scheme, with PME for
436 * Coulomb. It is not supported with only LJ PME. */
437 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
438 && ir->cutoff_scheme != ecutsGROUP);
440 pme_load_balancing_t* pme_loadbal = nullptr;
443 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata,
444 fr->nbv->useGpu(), &bPMETunePrinting);
447 if (!ir->bContinuation)
449 if (state->flags & (1U << estV))
451 auto v = makeArrayRef(state->v);
452 /* Set the velocities of vsites, shells and frozen atoms to zero */
453 for (i = 0; i < mdatoms->homenr; i++)
455 if (mdatoms->ptype[i] == eptVSite || mdatoms->ptype[i] == eptShell)
459 else if (mdatoms->cFREEZE)
461 for (m = 0; m < DIM; m++)
463 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
474 /* Constrain the initial coordinates and velocities */
475 do_constrain_first(fplog, constr, ir, mdatoms, state->natoms, state->x.arrayRefWithPadding(),
476 state->v.arrayRefWithPadding(), state->box, state->lambda[efptBONDED]);
480 /* Construct the virtual sites for the initial configuration */
481 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, nullptr, top.idef.iparams,
482 top.idef.il, fr->ePBC, fr->bMolPBC, cr, state->box);
486 if (ir->efep != efepNO)
488 /* Set free energy calculation frequency as the greatest common
489 * denominator of nstdhdl and repl_ex_nst. */
490 nstfep = ir->fepvals->nstdhdl;
493 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
495 if (useReplicaExchange)
497 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
501 /* Be REALLY careful about what flags you set here. You CANNOT assume
502 * this is the first step, since we might be restarting from a checkpoint,
503 * and in that case we should not do any modifications to the state.
505 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
507 // When restarting from a checkpoint, it can be appropriate to
508 // initialize ekind from quantities in the checkpoint. Otherwise,
509 // compute_globals must initialize ekind before the simulation
510 // starts/restarts. However, only the master rank knows what was
511 // found in the checkpoint file, so we have to communicate in
512 // order to coordinate the restart.
514 // TODO Consider removing this communication if/when checkpoint
515 // reading directly follows .tpr reading, because all ranks can
516 // agree on hasReadEkinState at that time.
517 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
520 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr);
522 if (hasReadEkinState)
524 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
527 unsigned int cglo_flags =
528 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
529 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
531 bSumEkinhOld = FALSE;
533 t_vcm vcm(top_global->groups, *ir);
534 reportComRemovalInfo(fplog, vcm);
536 /* To minimize communication, compute_globals computes the COM velocity
537 * and the kinetic energy for the velocities without COM motion removed.
538 * Thus to get the kinetic energy without the COM contribution, we need
539 * to call compute_globals twice.
541 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
543 unsigned int cglo_flags_iteration = cglo_flags;
544 if (bStopCM && cgloIteration == 0)
546 cglo_flags_iteration |= CGLO_STOPCM;
547 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
549 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
550 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, nullptr, enerd,
551 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
552 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
554 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
556 if (cglo_flags_iteration & CGLO_STOPCM)
558 /* At initialization, do not pass x with acceleration-correction mode
559 * to avoid (incorrect) correction of the initial coordinates.
561 rvec* xPtr = nullptr;
562 if (vcm.mode != ecmLINEAR_ACCELERATION_CORRECTION)
564 xPtr = state->x.rvec_array();
566 process_and_stopcm_grp(fplog, &vcm, *mdatoms, xPtr, state->v.rvec_array());
567 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
570 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
571 state->x.rvec_array(), state->box,
572 &shouldCheckNumberOfBondedInteractions);
573 if (ir->eI == eiVVAK)
575 /* a second call to get the half step temperature initialized as well */
576 /* we do the same call as above, but turn the pressure off -- internally to
577 compute_globals, this is recognized as a velocity verlet half-step
578 kinetic energy calculation. This minimized excess variables, but
579 perhaps loses some logic?*/
581 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
582 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, nullptr, enerd,
583 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
584 state->box, nullptr, &bSumEkinhOld, cglo_flags & ~CGLO_PRESSURE);
587 /* Calculate the initial half step temperature, and save the ekinh_old */
588 if (startingBehavior == StartingBehavior::NewSimulation)
590 for (i = 0; (i < ir->opts.ngtc); i++)
592 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
596 /* need to make an initiation call to get the Trotter variables set, as well as other constants
597 for non-trotter temperature control */
598 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
602 if (!ir->bContinuation)
604 if (constr && ir->eConstrAlg == econtLINCS)
606 fprintf(fplog, "RMS relative constraint deviation after constraining: %.2e\n",
609 if (EI_STATE_VELOCITY(ir->eI))
611 real temp = enerd->term[F_TEMP];
614 /* Result of Ekin averaged over velocities of -half
615 * and +half step, while we only have -half step here.
619 fprintf(fplog, "Initial temperature: %g K\n", temp);
624 fprintf(stderr, "starting mdrun '%s'\n", *(top_global->name));
627 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
631 sprintf(tbuf, "%s", "infinite");
633 if (ir->init_step > 0)
635 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
636 gmx_step_str(ir->init_step + ir->nsteps, sbuf), tbuf,
637 gmx_step_str(ir->init_step, sbuf2), ir->init_step * ir->delta_t);
641 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
643 fprintf(fplog, "\n");
646 walltime_accounting_start_time(walltime_accounting);
647 wallcycle_start(wcycle, ewcRUN);
648 print_start(fplog, cr, walltime_accounting, "mdrun");
651 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
652 int chkpt_ret = fcCheckPointParallel(cr->nodeid, NULL, 0);
655 gmx_fatal(3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0);
659 /***********************************************************
663 ************************************************************/
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;
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 =
678 (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
679 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
681 // Replica exchange, ensemble restraints and AWH need all
682 // simulations to remain synchronized, so they need
683 // checkpoints and stop conditions to act on the same step, so
684 // the propagation of such signals must take place between
685 // simulations, not just within simulations.
686 // TODO: Make algorithm initializers set these flags.
687 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
689 if (simulationsShareState)
691 // Inter-simulation signal communication does not need to happen
692 // often, so we use a minimum of 200 steps to reduce overhead.
693 const int c_minimumInterSimulationSignallingInterval = 200;
694 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
699 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
700 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
701 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
702 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
704 auto checkpointHandler = std::make_unique<CheckpointHandler>(
705 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]), simulationsShareState,
706 ir->nstlist == 0, MASTER(cr), mdrunOptions.writeConfout,
707 mdrunOptions.checkpointOptions.period);
709 const bool resetCountersIsLocal = true;
710 auto resetHandler = std::make_unique<ResetHandler>(
711 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
712 !resetCountersIsLocal, ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
713 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
715 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
717 step = ir->init_step;
720 // TODO extract this to new multi-simulation module
721 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
723 if (!multisim_int_all_are_equal(ms, ir->nsteps))
725 GMX_LOG(mdlog.warning)
727 "Note: The number of steps is not consistent across multi "
729 "but we are proceeding anyway!");
731 if (!multisim_int_all_are_equal(ms, ir->init_step))
733 GMX_LOG(mdlog.warning)
735 "Note: The initial step is not consistent across multi simulations,\n"
736 "but we are proceeding anyway!");
740 /* and stop now if we should */
741 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
745 /* Determine if this is a neighbor search step */
746 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
748 if (bPMETune && bNStList)
750 // This has to be here because PME load balancing is called so early.
751 // TODO: Move to after all booleans are defined.
752 if (useGpuForUpdate && !bFirstStep)
754 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
755 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
757 /* PME grid + cut-off optimization with GPUs or PME nodes */
758 pme_loadbal_do(pme_loadbal, cr, (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
759 fplog, mdlog, *ir, fr, state->box, state->x, wcycle, step, step_rel,
763 wallcycle_start(wcycle, ewcSTEP);
765 bLastStep = (step_rel == ir->nsteps);
766 t = t0 + step * ir->delta_t;
768 // TODO Refactor this, so that nstfep does not need a default value of zero
769 if (ir->efep != efepNO || ir->bSimTemp)
771 /* find and set the current lambdas */
772 setCurrentLambdasLocal(step, ir->fepvals, lam0, state->lambda, state->fep_state);
774 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
775 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
776 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
777 && (step > 0) && (startingBehavior == StartingBehavior::NewSimulation));
780 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
781 && do_per_step(step, replExParams.exchangeInterval));
783 if (doSimulatedAnnealing)
785 update_annealing_target_temp(ir, t, &upd);
788 /* Stop Center of Mass motion */
789 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
791 /* Determine whether or not to do Neighbour Searching */
792 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
794 /* Note that the stopHandler will cause termination at nstglobalcomm
795 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
796 * nstpcouple steps, we have computed the half-step kinetic energy
797 * of the previous step and can always output energies at the last step.
799 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
801 /* do_log triggers energy and virial calculation. Because this leads
802 * to different code paths, forces can be different. Thus for exact
803 * continuation we should avoid extra log output.
804 * Note that the || bLastStep can result in non-exact continuation
805 * beyond the last step. But we don't consider that to be an issue.
807 do_log = (do_per_step(step, ir->nstlog)
808 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
809 do_verbose = mdrunOptions.verbose
810 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
812 if (useGpuForUpdate && !bFirstStep && bNS)
814 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
815 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
816 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
817 // Copy coordinate from the GPU when needed at the search step.
818 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
819 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
820 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
821 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
824 if (bNS && !(bFirstStep && ir->bContinuation))
826 bMasterState = FALSE;
827 /* Correct the new box if it is too skewed */
828 if (inputrecDynamicBox(ir))
830 if (correct_box(fplog, step, state->box, graph))
833 // If update is offloaded, it should be informed about the box size change
837 set_pbc(&pbc, epbcXYZ, state->box);
838 integrator->setPbc(&pbc);
842 if (DOMAINDECOMP(cr) && bMasterState)
844 dd_collect_state(cr->dd, state, state_global);
847 if (DOMAINDECOMP(cr))
849 /* Repartition the domain decomposition */
850 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
851 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
852 fr, vsite, constr, nrnb, wcycle, do_verbose && !bPMETunePrinting);
853 shouldCheckNumberOfBondedInteractions = true;
854 upd.setNumAtoms(state->natoms);
858 if (MASTER(cr) && do_log)
860 energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
863 if (ir->efep != efepNO)
865 update_mdatoms(mdatoms, state->lambda[efptMASS]);
871 /* We need the kinetic energy at minus the half step for determining
872 * the full step kinetic energy and possibly for T-coupling.*/
873 /* This may not be quite working correctly yet . . . . */
874 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
875 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
876 nullptr, nullptr, nullptr, nullptr, mu_tot, constr, &nullSignaller,
877 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
878 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
879 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
880 &top, state->x.rvec_array(), state->box,
881 &shouldCheckNumberOfBondedInteractions);
883 clear_mat(force_vir);
885 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
887 /* Determine the energy and pressure:
888 * at nstcalcenergy steps and at energy output steps (set below).
890 if (EI_VV(ir->eI) && (!bInitStep))
892 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
893 bCalcVir = bCalcEnerStep
895 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
899 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
900 bCalcVir = bCalcEnerStep || (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
902 bCalcEner = bCalcEnerStep;
904 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
906 if (do_ene || do_log || bDoReplEx)
912 /* Do we need global communication ? */
913 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
914 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
916 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
917 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
918 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
922 /* Now is the time to relax the shells */
923 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
924 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
925 state->natoms, state->x.arrayRefWithPadding(),
926 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
927 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, graph,
928 shellfc, fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
932 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
933 is updated (or the AWH update will be performed twice for one step when continuing).
934 It would be best to call this update function from do_md_trajectory_writing but that
935 would occur after do_force. One would have to divide the update_awh function into one
936 function applying the AWH force and one doing the AWH bias update. The update AWH
937 bias function could then be called after do_md_trajectory_writing (then containing
938 update_awh_history). The checkpointing will in the future probably moved to the start
939 of the md loop which will rid of this issue. */
940 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
942 awh->updateHistory(state_global->awhHistory.get());
945 /* The coordinates (x) are shifted (to get whole molecules)
947 * This is parallellized as well, and does communication too.
948 * Check comments in sim_util.c
950 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation, imdSession, pull_work, step,
951 nrnb, wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
952 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, graph,
953 fr, runScheduleWork, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
954 (bNS ? GMX_FORCE_NS : 0) | force_flags, ddBalanceRegionHandler);
957 // VV integrators do not need the following velocity half step
958 // if it is the first step after starting from a checkpoint.
959 // That is, the half step is needed on all other steps, and
960 // also the first step when starting from a .tpr file.
961 if (EI_VV(ir->eI) && (!bFirstStep || startingBehavior == StartingBehavior::NewSimulation))
962 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
964 rvec* vbuf = nullptr;
966 wallcycle_start(wcycle, ewcUPDATE);
967 if (ir->eI == eiVV && bInitStep)
969 /* if using velocity verlet with full time step Ekin,
970 * take the first half step only to compute the
971 * virial for the first step. From there,
972 * revert back to the initial coordinates
973 * so that the input is actually the initial step.
975 snew(vbuf, state->natoms);
976 copy_rvecn(state->v.rvec_array(), vbuf, 0,
977 state->natoms); /* should make this better for parallelizing? */
981 /* this is for NHC in the Ekin(t+dt/2) version of vv */
982 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
983 trotter_seq, ettTSEQ1);
986 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
987 etrtVELOCITY1, cr, constr);
989 wallcycle_stop(wcycle, ewcUPDATE);
990 constrain_velocities(step, nullptr, state, shake_vir, constr, bCalcVir, do_log, do_ene);
991 wallcycle_start(wcycle, ewcUPDATE);
992 /* if VV, compute the pressure and constraints */
993 /* For VV2, we strictly only need this if using pressure
994 * control, but we really would like to have accurate pressures
996 * Think about ways around this in the future?
997 * For now, keep this choice in comments.
999 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1000 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1002 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1003 if (bCalcEner && ir->eI == eiVVAK)
1005 bSumEkinhOld = TRUE;
1007 /* for vv, the first half of the integration actually corresponds to the previous step.
1008 So we need information from the last step in the first half of the integration */
1009 if (bGStat || do_per_step(step - 1, nstglobalcomm))
1011 wallcycle_stop(wcycle, ewcUPDATE);
1012 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1013 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1014 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
1015 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1016 (bGStat ? CGLO_GSTAT : 0) | (bCalcEner ? CGLO_ENERGY : 0)
1017 | (bTemp ? CGLO_TEMPERATURE : 0) | (bPres ? CGLO_PRESSURE : 0)
1018 | (bPres ? CGLO_CONSTRAINT : 0) | (bStopCM ? CGLO_STOPCM : 0)
1019 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
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);
1034 process_and_stopcm_grp(fplog, &vcm, *mdatoms, state->x.rvec_array(),
1035 state->v.rvec_array());
1036 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1038 wallcycle_start(wcycle, ewcUPDATE);
1040 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1045 m_add(force_vir, shake_vir,
1046 total_vir); /* we need the un-dispersion corrected total vir here */
1047 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1048 trotter_seq, ettTSEQ2);
1050 /* TODO This is only needed when we're about to write
1051 * a checkpoint, because we use it after the restart
1052 * (in a kludge?). But what should we be doing if
1053 * the startingBehavior is NewSimulation or bInitStep are true? */
1054 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1056 copy_mat(shake_vir, state->svir_prev);
1057 copy_mat(force_vir, state->fvir_prev);
1059 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1061 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1062 enerd->term[F_TEMP] =
1063 sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1064 enerd->term[F_EKIN] = trace(ekind->ekin);
1067 else if (bExchanged)
1069 wallcycle_stop(wcycle, ewcUPDATE);
1070 /* We need the kinetic energy at minus the half step for determining
1071 * the full step kinetic energy and possibly for T-coupling.*/
1072 /* This may not be quite working correctly yet . . . . */
1073 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(),
1074 state->v.rvec_array(), state->box, state->lambda[efptVDW],
1075 mdatoms, nrnb, &vcm, wcycle, enerd, nullptr, nullptr, nullptr,
1076 nullptr, mu_tot, constr, &nullSignaller, state->box, nullptr,
1077 &bSumEkinhOld, CGLO_GSTAT | CGLO_TEMPERATURE);
1078 wallcycle_start(wcycle, ewcUPDATE);
1081 /* if it's the initial step, we performed this first step just to get the constraint virial */
1082 if (ir->eI == eiVV && bInitStep)
1084 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1087 wallcycle_stop(wcycle, ewcUPDATE);
1090 /* compute the conserved quantity */
1093 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1096 last_ekin = enerd->term[F_EKIN];
1098 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1100 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1102 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1103 if (ir->efep != efepNO)
1105 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1109 /* ######## END FIRST UPDATE STEP ############## */
1110 /* ######## If doing VV, we now have v(dt) ###### */
1113 /* perform extended ensemble sampling in lambda - we don't
1114 actually move to the new state before outputting
1115 statistics, but if performing simulated tempering, we
1116 do update the velocities and the tau_t. */
1118 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state,
1119 state->dfhist, step, state->v.rvec_array(), mdatoms);
1120 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1123 copy_df_history(state_global->dfhist, state->dfhist);
1127 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1128 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1129 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1130 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1131 || checkpointHandler->isCheckpointingStep()))
1133 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1134 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1136 // Copy velocities if needed for the output/checkpointing.
1137 // NOTE: Copy on the search steps is done at the beginning of the step.
1138 if (useGpuForUpdate && !bNS
1139 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1141 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1142 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1144 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1145 // and update is offloaded hence forces are kept on the GPU for update and have not been
1146 // already transferred in do_force().
1147 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1148 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1149 // prior to GPU update.
1150 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1151 // copy call in do_force(...).
1152 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1153 // on host after the D2H copy in do_force(...).
1154 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1155 && do_per_step(step, ir->nstfout))
1157 stateGpu->copyForcesFromGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1158 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1160 /* Now we have the energies and forces corresponding to the
1161 * coordinates at time t. We must output all of this before
1164 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state, state_global,
1165 observablesHistory, top_global, fr, outf, energyOutput, ekind, f,
1166 checkpointHandler->isCheckpointingStep(), bRerunMD, bLastStep,
1167 mdrunOptions.writeConfout, bSumEkinhOld);
1168 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1169 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1171 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1172 if (startingBehavior != StartingBehavior::NewSimulation
1173 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1175 copy_mat(state->svir_prev, shake_vir);
1176 copy_mat(state->fvir_prev, force_vir);
1179 stopHandler->setSignal();
1180 resetHandler->setSignal(walltime_accounting);
1182 if (bGStat || !PAR(cr))
1184 /* In parallel we only have to check for checkpointing in steps
1185 * where we do global communication,
1186 * otherwise the other nodes don't know.
1188 checkpointHandler->setSignal(walltime_accounting);
1191 /* ######### START SECOND UPDATE STEP ################# */
1193 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1194 controlled in preprocessing */
1196 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1198 gmx_bool bIfRandomize;
1199 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1200 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1201 if (constr && bIfRandomize)
1203 constrain_velocities(step, nullptr, state, tmp_vir, constr, bCalcVir, do_log, do_ene);
1206 /* Box is changed in update() when we do pressure coupling,
1207 * but we should still use the old box for energy corrections and when
1208 * writing it to the energy file, so it matches the trajectory files for
1209 * the same timestep above. Make a copy in a separate array.
1211 copy_mat(state->box, lastbox);
1215 wallcycle_start(wcycle, ewcUPDATE);
1216 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1219 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1220 /* We can only do Berendsen coupling after we have summed
1221 * the kinetic energy or virial. Since the happens
1222 * in global_state after update, we should only do it at
1223 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1228 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1229 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1234 /* velocity half-step update */
1235 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1236 etrtVELOCITY2, cr, constr);
1239 /* Above, initialize just copies ekinh into ekin,
1240 * it doesn't copy position (for VV),
1241 * and entire integrator for MD.
1244 if (ir->eI == eiVVAK)
1246 cbuf.resize(state->x.size());
1247 std::copy(state->x.begin(), state->x.end(), cbuf.begin());
1250 /* With leap-frog type integrators we compute the kinetic energy
1251 * at a whole time step as the average of the half-time step kinetic
1252 * energies of two subsequent steps. Therefore we need to compute the
1253 * half step kinetic energy also if we need energies at the next step.
1255 const bool needHalfStepKineticEnergy =
1256 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1258 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1259 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1260 const bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN
1261 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1263 if (useGpuForUpdate)
1267 integrator->set(stateGpu->getCoordinates(), stateGpu->getVelocities(),
1268 stateGpu->getForces(), top.idef, *mdatoms, ekind->ngtc);
1270 // Copy data to the GPU after buffers might have being reinitialized
1271 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1272 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1275 // If the buffer ops were not offloaded this step, the forces are on the host and have to be copied
1276 if (!runScheduleWork->stepWork.useGpuFBufferOps)
1278 stateGpu->copyForcesToGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1281 const bool doTemperatureScaling =
1282 (ir->etc != etcNO && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1284 // This applies Leap-Frog, LINCS and SETTLE in succession
1285 integrator->integrate(stateGpu->getForcesReadyOnDeviceEvent(
1286 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1287 ir->delta_t, true, bCalcVir, shake_vir, doTemperatureScaling,
1288 ekind->tcstat, doParrinelloRahman, ir->nstpcouple * ir->delta_t, M);
1290 // Copy velocities D2H after update if:
1291 // - Globals are computed this step (includes the energy output steps).
1292 // - Temperature is needed for the next step.
1293 if (bGStat || needHalfStepKineticEnergy)
1295 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1296 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1301 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1302 etrtPOSITION, cr, constr);
1304 wallcycle_stop(wcycle, ewcUPDATE);
1306 constrain_coordinates(step, &dvdl_constr, state, shake_vir, &upd, constr, bCalcVir,
1309 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state, cr, nrnb, wcycle, &upd,
1310 constr, do_log, do_ene);
1311 finish_update(ir, mdatoms, state, graph, nrnb, wcycle, &upd, constr);
1314 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1316 updatePrevStepPullCom(pull_work, state);
1319 if (ir->eI == eiVVAK)
1321 /* erase F_EKIN and F_TEMP here? */
1322 /* just compute the kinetic energy at the half step to perform a trotter step */
1323 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1324 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1325 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller, lastbox,
1326 nullptr, &bSumEkinhOld, (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE);
1327 wallcycle_start(wcycle, ewcUPDATE);
1328 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1329 /* now we know the scaling, we can compute the positions again */
1330 std::copy(cbuf.begin(), cbuf.end(), state->x.begin());
1332 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1333 etrtPOSITION, cr, constr);
1334 wallcycle_stop(wcycle, ewcUPDATE);
1336 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1337 /* are the small terms in the shake_vir here due
1338 * to numerical errors, or are they important
1339 * physically? I'm thinking they are just errors, but not completely sure.
1340 * For now, will call without actually constraining, constr=NULL*/
1341 finish_update(ir, mdatoms, state, graph, nrnb, wcycle, &upd, nullptr);
1345 /* this factor or 2 correction is necessary
1346 because half of the constraint force is removed
1347 in the vv step, so we have to double it. See
1348 the Redmine issue #1255. It is not yet clear
1349 if the factor of 2 is exact, or just a very
1350 good approximation, and this will be
1351 investigated. The next step is to see if this
1352 can be done adding a dhdl contribution from the
1353 rattle step, but this is somewhat more
1354 complicated with the current code. Will be
1355 investigated, hopefully for 4.6.3. However,
1356 this current solution is much better than
1357 having it completely wrong.
1359 enerd->term[F_DVDL_CONSTR] += 2 * dvdl_constr;
1363 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1366 if (vsite != nullptr)
1368 wallcycle_start(wcycle, ewcVSITECONSTR);
1369 if (graph != nullptr)
1371 shift_self(graph, state->box, state->x.rvec_array());
1373 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, state->v.rvec_array(),
1374 top.idef.iparams, top.idef.il, fr->ePBC, fr->bMolPBC, cr, state->box);
1376 if (graph != nullptr)
1378 unshift_self(graph, state->box, state->x.rvec_array());
1380 wallcycle_stop(wcycle, ewcVSITECONSTR);
1383 /* ############## IF NOT VV, Calculate globals HERE ############ */
1384 /* With Leap-Frog we can skip compute_globals at
1385 * non-communication steps, but we need to calculate
1386 * the kinetic energy one step before communication.
1389 // Organize to do inter-simulation signalling on steps if
1390 // and when algorithms require it.
1391 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1393 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1395 // Copy coordinates when needed to stop the CM motion.
1396 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1398 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1399 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1401 // Since we're already communicating at this step, we
1402 // can propagate intra-simulation signals. Note that
1403 // check_nstglobalcomm has the responsibility for
1404 // choosing the value of nstglobalcomm that is one way
1405 // bGStat becomes true, so we can't get into a
1406 // situation where e.g. checkpointing can't be
1408 bool doIntraSimSignal = true;
1409 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1412 gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1413 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1414 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &signaller, lastbox,
1415 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1416 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1417 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1418 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1419 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1420 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1422 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1423 top_global, &top, state->x.rvec_array(), state->box,
1424 &shouldCheckNumberOfBondedInteractions);
1425 if (!EI_VV(ir->eI) && bStopCM)
1427 process_and_stopcm_grp(fplog, &vcm, *mdatoms, state->x.rvec_array(),
1428 state->v.rvec_array());
1429 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1431 // TODO: The special case of removing CM motion should be dealt more gracefully
1432 if (useGpuForUpdate)
1434 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1435 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1441 /* ############# END CALC EKIN AND PRESSURE ################# */
1443 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1444 the virial that should probably be addressed eventually. state->veta has better properies,
1445 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1446 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1448 if (ir->efep != efepNO && !EI_VV(ir->eI))
1450 /* Sum up the foreign energy and dhdl terms for md and sd.
1451 Currently done every step so that dhdl is correct in the .edr */
1452 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1455 update_pcouple_after_coordinates(fplog, step, ir, mdatoms, pres, force_vir, shake_vir,
1456 pressureCouplingMu, state, nrnb, &upd, !useGpuForUpdate);
1458 const bool doBerendsenPressureCoupling =
1459 (inputrec->epc == epcBERENDSEN && do_per_step(step, inputrec->nstpcouple));
1460 if (useGpuForUpdate && (doBerendsenPressureCoupling || doParrinelloRahman))
1462 integrator->scaleCoordinates(pressureCouplingMu);
1464 set_pbc(&pbc, epbcXYZ, state->box);
1465 integrator->setPbc(&pbc);
1468 /* ################# END UPDATE STEP 2 ################# */
1469 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1471 /* The coordinates (x) were unshifted in update */
1474 /* We will not sum ekinh_old,
1475 * so signal that we still have to do it.
1477 bSumEkinhOld = TRUE;
1482 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1484 /* use the directly determined last velocity, not actually the averaged half steps */
1485 if (bTrotter && ir->eI == eiVV)
1487 enerd->term[F_EKIN] = last_ekin;
1489 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1491 if (integratorHasConservedEnergyQuantity(ir))
1495 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1499 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1502 /* ######### END PREPARING EDR OUTPUT ########### */
1508 if (fplog && do_log && bDoExpanded)
1510 /* only needed if doing expanded ensemble */
1511 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals,
1512 ir->bSimTemp ? ir->simtempvals : nullptr,
1513 state_global->dfhist, state->fep_state, ir->nstlog, step);
1517 energyOutput.addDataAtEnergyStep(bDoDHDL, bCalcEnerStep, t, mdatoms->tmass, enerd, state,
1518 ir->fepvals, ir->expandedvals, lastbox, shake_vir,
1519 force_vir, total_vir, pres, ekind, mu_tot, constr);
1523 energyOutput.recordNonEnergyStep();
1526 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1527 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1529 if (doSimulatedAnnealing)
1531 energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
1533 if (do_log || do_ene || do_dr || do_or)
1535 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1536 do_log ? fplog : nullptr, step, t, fcd, awh.get());
1541 pull_print_output(pull_work, step, t);
1544 if (do_per_step(step, ir->nstlog))
1546 if (fflush(fplog) != 0)
1548 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1554 /* Have to do this part _after_ outputting the logfile and the edr file */
1555 /* Gets written into the state at the beginning of next loop*/
1556 state->fep_state = lamnew;
1558 /* Print the remaining wall clock time for the run */
1559 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1563 fprintf(stderr, "\n");
1565 print_time(stderr, walltime_accounting, step, ir, cr);
1568 /* Ion/water position swapping.
1569 * Not done in last step since trajectory writing happens before this call
1570 * in the MD loop and exchanges would be lost anyway. */
1571 bNeedRepartition = FALSE;
1572 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep && do_per_step(step, ir->swap->nstswap))
1575 do_swapcoords(cr, step, t, ir, swap, wcycle, as_rvec_array(state->x.data()),
1576 state->box, MASTER(cr) && mdrunOptions.verbose, bRerunMD);
1578 if (bNeedRepartition && DOMAINDECOMP(cr))
1580 dd_collect_state(cr->dd, state, state_global);
1584 /* Replica exchange */
1588 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1591 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1593 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1, state_global, *top_global, ir,
1594 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
1595 nrnb, wcycle, FALSE);
1596 shouldCheckNumberOfBondedInteractions = true;
1597 upd.setNumAtoms(state->natoms);
1603 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1604 /* With all integrators, except VV, we need to retain the pressure
1605 * at the current step for coupling at the next step.
1607 if ((state->flags & (1U << estPRES_PREV))
1608 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1610 /* Store the pressure in t_state for pressure coupling
1611 * at the next MD step.
1613 copy_mat(pres, state->pres_prev);
1616 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1618 if ((membed != nullptr) && (!bLastStep))
1620 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1623 cycles = wallcycle_stop(wcycle, ewcSTEP);
1624 if (DOMAINDECOMP(cr) && wcycle)
1626 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1629 /* increase the MD step number */
1633 resetHandler->resetCounters(step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb,
1634 fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1636 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1637 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1639 /* End of main MD loop */
1641 /* Closing TNG files can include compressing data. Therefore it is good to do that
1642 * before stopping the time measurements. */
1643 mdoutf_tng_close(outf);
1645 /* Stop measuring walltime */
1646 walltime_accounting_end_time(walltime_accounting);
1648 if (!thisRankHasDuty(cr, DUTY_PME))
1650 /* Tell the PME only node to finish */
1651 gmx_pme_send_finish(cr);
1656 if (ir->nstcalcenergy > 0)
1658 energyOutput.printAnnealingTemperatures(fplog, groups, &(ir->opts));
1659 energyOutput.printAverages(fplog, groups);
1666 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1669 done_shellfc(fplog, shellfc, step_rel);
1671 if (useReplicaExchange && MASTER(cr))
1673 print_replica_exchange_statistics(fplog, repl_ex);
1676 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1678 global_stat_destroy(gstat);