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)
338 || (simulationWork.useGpuDirectCommunication
339 && simulationWork.useGpuPmePpCommunication),
340 "Domain decomposition is not supported with the GPU update when not "
341 "using direct GPU communication.\n");
342 GMX_RELEASE_ASSERT(useGpuForPme || (useGpuForNonbonded && simulationWork.useGpuBufferOps),
343 "Either PME or short-ranged non-bonded interaction tasks must run on "
344 "the GPU to use GPU update.\n");
345 GMX_RELEASE_ASSERT(ir->eI == eiMD,
346 "Only the md integrator is supported with the GPU update.\n");
348 ir->etc != etcNOSEHOOVER,
349 "Nose-Hoover temperature coupling is not supported with the GPU update.\n");
350 GMX_RELEASE_ASSERT(ir->epc == epcNO || ir->epc == epcPARRINELLORAHMAN || ir->epc == epcBERENDSEN,
351 "Only Parrinello-Rahman and Berendsen pressure coupling are supported "
352 "with the GPU update.\n");
353 GMX_RELEASE_ASSERT(!mdatoms->haveVsites,
354 "Virtual sites are not supported with the GPU update.\n");
355 GMX_RELEASE_ASSERT(ed == nullptr,
356 "Essential dynamics is not supported with the GPU update.\n");
357 GMX_RELEASE_ASSERT(!ir->bPull && !ir->pull,
358 "Pulling is not supported with the GPU update.\n");
359 GMX_RELEASE_ASSERT(fcd->orires.nr == 0,
360 "Orientation restraints are not supported with the GPU update.\n");
361 GMX_RELEASE_ASSERT(ir->efep == efepNO,
362 "Free energy perturbations are not supported with the GPU update.");
363 GMX_RELEASE_ASSERT(graph == nullptr, "The graph is not supported with GPU update.");
365 if (constr != nullptr && constr->numConstraintsTotal() > 0)
369 .appendText("Updating coordinates and applying constraints on the GPU.");
373 GMX_LOG(mdlog.info).asParagraph().appendText("Updating coordinates on the GPU.");
375 integrator = std::make_unique<UpdateConstrainCuda>(
376 *ir, *top_global, stateGpu->getUpdateStream(), stateGpu->xUpdatedOnDevice());
379 set_pbc(&pbc, epbcXYZ, state->box);
380 integrator->setPbc(&pbc);
383 if (useGpuForPme || (useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
385 changePinningPolicy(&state->x, PinningPolicy::PinnedIfSupported);
387 if ((useGpuForNonbonded && useGpuForBufferOps) || useGpuForUpdate)
389 changePinningPolicy(&f, PinningPolicy::PinnedIfSupported);
393 changePinningPolicy(&state->v, PinningPolicy::PinnedIfSupported);
395 /*****************************************************************************************/
397 // NOTE: The global state is no longer used at this point.
398 // But state_global is still used as temporary storage space for writing
399 // the global state to file and potentially for replica exchange.
400 // (Global topology should persist.)
402 update_mdatoms(mdatoms, state->lambda[efptMASS]);
406 /* Check nstexpanded here, because the grompp check was broken */
407 if (ir->expandedvals->nstexpanded % ir->nstcalcenergy != 0)
410 "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
412 init_expanded_ensemble(startingBehavior != StartingBehavior::NewSimulation, ir, state->dfhist);
417 EnergyElement::initializeEnergyHistory(startingBehavior, observablesHistory, &energyOutput);
420 preparePrevStepPullCom(ir, pull_work, mdatoms, state, state_global, cr,
421 startingBehavior != StartingBehavior::NewSimulation);
423 // TODO: Remove this by converting AWH into a ForceProvider
424 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms,
425 startingBehavior != StartingBehavior::NewSimulation,
426 shellfc != nullptr, opt2fn("-awh", nfile, fnm), pull_work);
428 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
429 if (useReplicaExchange && MASTER(cr))
431 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir, replExParams);
433 /* PME tuning is only supported in the Verlet scheme, with PME for
434 * Coulomb. It is not supported with only LJ PME. */
435 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
436 && ir->cutoff_scheme != ecutsGROUP);
438 pme_load_balancing_t* pme_loadbal = nullptr;
441 pme_loadbal_init(&pme_loadbal, cr, mdlog, *ir, state->box, *fr->ic, *fr->nbv, fr->pmedata,
442 fr->nbv->useGpu(), &bPMETunePrinting);
445 if (!ir->bContinuation)
447 if (state->flags & (1U << estV))
449 auto v = makeArrayRef(state->v);
450 /* Set the velocities of vsites, shells and frozen atoms to zero */
451 for (i = 0; i < mdatoms->homenr; i++)
453 if (mdatoms->ptype[i] == eptVSite || mdatoms->ptype[i] == eptShell)
457 else if (mdatoms->cFREEZE)
459 for (m = 0; m < DIM; m++)
461 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
472 /* Constrain the initial coordinates and velocities */
473 do_constrain_first(fplog, constr, ir, mdatoms, state->natoms, state->x.arrayRefWithPadding(),
474 state->v.arrayRefWithPadding(), state->box, state->lambda[efptBONDED]);
478 /* Construct the virtual sites for the initial configuration */
479 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, nullptr, top.idef.iparams,
480 top.idef.il, fr->ePBC, fr->bMolPBC, cr, state->box);
484 if (ir->efep != efepNO)
486 /* Set free energy calculation frequency as the greatest common
487 * denominator of nstdhdl and repl_ex_nst. */
488 nstfep = ir->fepvals->nstdhdl;
491 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
493 if (useReplicaExchange)
495 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
499 /* Be REALLY careful about what flags you set here. You CANNOT assume
500 * this is the first step, since we might be restarting from a checkpoint,
501 * and in that case we should not do any modifications to the state.
503 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
505 // When restarting from a checkpoint, it can be appropriate to
506 // initialize ekind from quantities in the checkpoint. Otherwise,
507 // compute_globals must initialize ekind before the simulation
508 // starts/restarts. However, only the master rank knows what was
509 // found in the checkpoint file, so we have to communicate in
510 // order to coordinate the restart.
512 // TODO Consider removing this communication if/when checkpoint
513 // reading directly follows .tpr reading, because all ranks can
514 // agree on hasReadEkinState at that time.
515 bool hasReadEkinState = MASTER(cr) ? state_global->ekinstate.hasReadEkinState : false;
518 gmx_bcast(sizeof(hasReadEkinState), &hasReadEkinState, cr);
520 if (hasReadEkinState)
522 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
525 unsigned int cglo_flags =
526 (CGLO_TEMPERATURE | CGLO_GSTAT | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
527 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0) | (hasReadEkinState ? CGLO_READEKIN : 0));
529 bSumEkinhOld = FALSE;
531 t_vcm vcm(top_global->groups, *ir);
532 reportComRemovalInfo(fplog, vcm);
534 /* To minimize communication, compute_globals computes the COM velocity
535 * and the kinetic energy for the velocities without COM motion removed.
536 * Thus to get the kinetic energy without the COM contribution, we need
537 * to call compute_globals twice.
539 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
541 unsigned int cglo_flags_iteration = cglo_flags;
542 if (bStopCM && cgloIteration == 0)
544 cglo_flags_iteration |= CGLO_STOPCM;
545 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
547 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
548 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, nullptr, enerd,
549 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
550 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
552 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
554 if (cglo_flags_iteration & CGLO_STOPCM)
556 /* At initialization, do not pass x with acceleration-correction mode
557 * to avoid (incorrect) correction of the initial coordinates.
559 rvec* xPtr = nullptr;
560 if (vcm.mode != ecmLINEAR_ACCELERATION_CORRECTION)
562 xPtr = state->x.rvec_array();
564 process_and_stopcm_grp(fplog, &vcm, *mdatoms, xPtr, state->v.rvec_array());
565 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
568 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
569 state->x.rvec_array(), state->box,
570 &shouldCheckNumberOfBondedInteractions);
571 if (ir->eI == eiVVAK)
573 /* a second call to get the half step temperature initialized as well */
574 /* we do the same call as above, but turn the pressure off -- internally to
575 compute_globals, this is recognized as a velocity verlet half-step
576 kinetic energy calculation. This minimized excess variables, but
577 perhaps loses some logic?*/
579 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
580 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, nullptr, enerd,
581 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
582 state->box, nullptr, &bSumEkinhOld, cglo_flags & ~CGLO_PRESSURE);
585 /* Calculate the initial half step temperature, and save the ekinh_old */
586 if (startingBehavior == StartingBehavior::NewSimulation)
588 for (i = 0; (i < ir->opts.ngtc); i++)
590 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
594 /* need to make an initiation call to get the Trotter variables set, as well as other constants
595 for non-trotter temperature control */
596 auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
600 if (!ir->bContinuation)
602 if (constr && ir->eConstrAlg == econtLINCS)
604 fprintf(fplog, "RMS relative constraint deviation after constraining: %.2e\n",
607 if (EI_STATE_VELOCITY(ir->eI))
609 real temp = enerd->term[F_TEMP];
612 /* Result of Ekin averaged over velocities of -half
613 * and +half step, while we only have -half step here.
617 fprintf(fplog, "Initial temperature: %g K\n", temp);
622 fprintf(stderr, "starting mdrun '%s'\n", *(top_global->name));
625 sprintf(tbuf, "%8.1f", (ir->init_step + ir->nsteps) * ir->delta_t);
629 sprintf(tbuf, "%s", "infinite");
631 if (ir->init_step > 0)
633 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
634 gmx_step_str(ir->init_step + ir->nsteps, sbuf), tbuf,
635 gmx_step_str(ir->init_step, sbuf2), ir->init_step * ir->delta_t);
639 fprintf(stderr, "%s steps, %s ps.\n", gmx_step_str(ir->nsteps, sbuf), tbuf);
641 fprintf(fplog, "\n");
644 walltime_accounting_start_time(walltime_accounting);
645 wallcycle_start(wcycle, ewcRUN);
646 print_start(fplog, cr, walltime_accounting, "mdrun");
649 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
650 int chkpt_ret = fcCheckPointParallel(cr->nodeid, NULL, 0);
653 gmx_fatal(3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0);
657 /***********************************************************
661 ************************************************************/
664 /* Skip the first Nose-Hoover integration when we get the state from tpx */
665 bInitStep = startingBehavior == StartingBehavior::NewSimulation || EI_VV(ir->eI);
666 bSumEkinhOld = FALSE;
668 bNeedRepartition = FALSE;
670 bool simulationsShareState = false;
671 int nstSignalComm = nstglobalcomm;
673 // TODO This implementation of ensemble orientation restraints is nasty because
674 // a user can't just do multi-sim with single-sim orientation restraints.
675 bool usingEnsembleRestraints =
676 (fcd->disres.nsystems > 1) || ((ms != nullptr) && (fcd->orires.nr != 0));
677 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && (ms != nullptr));
679 // Replica exchange, ensemble restraints and AWH need all
680 // simulations to remain synchronized, so they need
681 // checkpoints and stop conditions to act on the same step, so
682 // the propagation of such signals must take place between
683 // simulations, not just within simulations.
684 // TODO: Make algorithm initializers set these flags.
685 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
687 if (simulationsShareState)
689 // Inter-simulation signal communication does not need to happen
690 // often, so we use a minimum of 200 steps to reduce overhead.
691 const int c_minimumInterSimulationSignallingInterval = 200;
692 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1) / nstglobalcomm)
697 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
698 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), simulationsShareState,
699 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstSignalComm,
700 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
702 auto checkpointHandler = std::make_unique<CheckpointHandler>(
703 compat::make_not_null<SimulationSignal*>(&signals[eglsCHKPT]), simulationsShareState,
704 ir->nstlist == 0, MASTER(cr), mdrunOptions.writeConfout,
705 mdrunOptions.checkpointOptions.period);
707 const bool resetCountersIsLocal = true;
708 auto resetHandler = std::make_unique<ResetHandler>(
709 compat::make_not_null<SimulationSignal*>(&signals[eglsRESETCOUNTERS]),
710 !resetCountersIsLocal, ir->nsteps, MASTER(cr), mdrunOptions.timingOptions.resetHalfway,
711 mdrunOptions.maximumHoursToRun, mdlog, wcycle, walltime_accounting);
713 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
715 step = ir->init_step;
718 // TODO extract this to new multi-simulation module
719 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
721 if (!multisim_int_all_are_equal(ms, ir->nsteps))
723 GMX_LOG(mdlog.warning)
725 "Note: The number of steps is not consistent across multi "
727 "but we are proceeding anyway!");
729 if (!multisim_int_all_are_equal(ms, ir->init_step))
731 GMX_LOG(mdlog.warning)
733 "Note: The initial step is not consistent across multi simulations,\n"
734 "but we are proceeding anyway!");
738 /* and stop now if we should */
739 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
743 /* Determine if this is a neighbor search step */
744 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
746 if (bPMETune && bNStList)
748 // This has to be here because PME load balancing is called so early.
749 // TODO: Move to after all booleans are defined.
750 if (useGpuForUpdate && !bFirstStep)
752 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
753 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
755 /* PME grid + cut-off optimization with GPUs or PME nodes */
756 pme_loadbal_do(pme_loadbal, cr, (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
757 fplog, mdlog, *ir, fr, state->box, state->x, wcycle, step, step_rel,
761 wallcycle_start(wcycle, ewcSTEP);
763 bLastStep = (step_rel == ir->nsteps);
764 t = t0 + step * ir->delta_t;
766 // TODO Refactor this, so that nstfep does not need a default value of zero
767 if (ir->efep != efepNO || ir->bSimTemp)
769 /* find and set the current lambdas */
770 setCurrentLambdasLocal(step, ir->fepvals, lam0, state->lambda, state->fep_state);
772 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
773 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
774 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded) && (ir->bExpanded)
775 && (step > 0) && (startingBehavior == StartingBehavior::NewSimulation));
778 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep
779 && do_per_step(step, replExParams.exchangeInterval));
781 if (doSimulatedAnnealing)
783 update_annealing_target_temp(ir, t, &upd);
786 /* Stop Center of Mass motion */
787 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
789 /* Determine whether or not to do Neighbour Searching */
790 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
792 /* Note that the stopHandler will cause termination at nstglobalcomm
793 * steps. Since this concides with nstcalcenergy, nsttcouple and/or
794 * nstpcouple steps, we have computed the half-step kinetic energy
795 * of the previous step and can always output energies at the last step.
797 bLastStep = bLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
799 /* do_log triggers energy and virial calculation. Because this leads
800 * to different code paths, forces can be different. Thus for exact
801 * continuation we should avoid extra log output.
802 * Note that the || bLastStep can result in non-exact continuation
803 * beyond the last step. But we don't consider that to be an issue.
805 do_log = (do_per_step(step, ir->nstlog)
806 || (bFirstStep && startingBehavior == StartingBehavior::NewSimulation) || bLastStep);
807 do_verbose = mdrunOptions.verbose
808 && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
810 if (useGpuForUpdate && !bFirstStep && bNS)
812 // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
813 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
814 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
815 // Copy coordinate from the GPU when needed at the search step.
816 // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
817 // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
818 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
819 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
822 if (bNS && !(bFirstStep && ir->bContinuation))
824 bMasterState = FALSE;
825 /* Correct the new box if it is too skewed */
826 if (inputrecDynamicBox(ir))
828 if (correct_box(fplog, step, state->box, graph))
831 // If update is offloaded, it should be informed about the box size change
835 set_pbc(&pbc, epbcXYZ, state->box);
836 integrator->setPbc(&pbc);
840 if (DOMAINDECOMP(cr) && bMasterState)
842 dd_collect_state(cr->dd, state, state_global);
845 if (DOMAINDECOMP(cr))
847 /* Repartition the domain decomposition */
848 dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
849 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
850 fr, vsite, constr, nrnb, wcycle, do_verbose && !bPMETunePrinting);
851 shouldCheckNumberOfBondedInteractions = true;
852 upd.setNumAtoms(state->natoms);
856 if (MASTER(cr) && do_log)
858 energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
861 if (ir->efep != efepNO)
863 update_mdatoms(mdatoms, state->lambda[efptMASS]);
869 /* We need the kinetic energy at minus the half step for determining
870 * the full step kinetic energy and possibly for T-coupling.*/
871 /* This may not be quite working correctly yet . . . . */
872 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
873 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
874 nullptr, nullptr, nullptr, nullptr, mu_tot, constr, &nullSignaller,
875 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
876 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
877 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
878 &top, state->x.rvec_array(), state->box,
879 &shouldCheckNumberOfBondedInteractions);
881 clear_mat(force_vir);
883 checkpointHandler->decideIfCheckpointingThisStep(bNS, bFirstStep, bLastStep);
885 /* Determine the energy and pressure:
886 * at nstcalcenergy steps and at energy output steps (set below).
888 if (EI_VV(ir->eI) && (!bInitStep))
890 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
891 bCalcVir = bCalcEnerStep
893 && (do_per_step(step, ir->nstpcouple) || do_per_step(step - 1, ir->nstpcouple)));
897 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
898 bCalcVir = bCalcEnerStep || (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
900 bCalcEner = bCalcEnerStep;
902 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
904 if (do_ene || do_log || bDoReplEx)
910 /* Do we need global communication ? */
911 bGStat = (bCalcVir || bCalcEner || bStopCM || do_per_step(step, nstglobalcomm)
912 || (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step - 1, nstglobalcomm)));
914 force_flags = (GMX_FORCE_STATECHANGED | ((inputrecDynamicBox(ir)) ? GMX_FORCE_DYNAMICBOX : 0)
915 | GMX_FORCE_ALLFORCES | (bCalcVir ? GMX_FORCE_VIRIAL : 0)
916 | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
920 /* Now is the time to relax the shells */
921 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
922 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
923 state->natoms, state->x.arrayRefWithPadding(),
924 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
925 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, graph,
926 shellfc, fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
930 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias
931 is updated (or the AWH update will be performed twice for one step when continuing).
932 It would be best to call this update function from do_md_trajectory_writing but that
933 would occur after do_force. One would have to divide the update_awh function into one
934 function applying the AWH force and one doing the AWH bias update. The update AWH
935 bias function could then be called after do_md_trajectory_writing (then containing
936 update_awh_history). The checkpointing will in the future probably moved to the start
937 of the md loop which will rid of this issue. */
938 if (awh && checkpointHandler->isCheckpointingStep() && MASTER(cr))
940 awh->updateHistory(state_global->awhHistory.get());
943 /* The coordinates (x) are shifted (to get whole molecules)
945 * This is parallellized as well, and does communication too.
946 * Check comments in sim_util.c
948 do_force(fplog, cr, ms, ir, awh.get(), enforcedRotation, imdSession, pull_work, step,
949 nrnb, wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
950 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, graph,
951 fr, runScheduleWork, vsite, mu_tot, t, ed ? ed->getLegacyED() : nullptr,
952 (bNS ? GMX_FORCE_NS : 0) | force_flags, ddBalanceRegionHandler);
955 // VV integrators do not need the following velocity half step
956 // if it is the first step after starting from a checkpoint.
957 // That is, the half step is needed on all other steps, and
958 // also the first step when starting from a .tpr file.
959 if (EI_VV(ir->eI) && (!bFirstStep || startingBehavior == StartingBehavior::NewSimulation))
960 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
962 rvec* vbuf = nullptr;
964 wallcycle_start(wcycle, ewcUPDATE);
965 if (ir->eI == eiVV && bInitStep)
967 /* if using velocity verlet with full time step Ekin,
968 * take the first half step only to compute the
969 * virial for the first step. From there,
970 * revert back to the initial coordinates
971 * so that the input is actually the initial step.
973 snew(vbuf, state->natoms);
974 copy_rvecn(state->v.rvec_array(), vbuf, 0,
975 state->natoms); /* should make this better for parallelizing? */
979 /* this is for NHC in the Ekin(t+dt/2) version of vv */
980 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
981 trotter_seq, ettTSEQ1);
984 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
985 etrtVELOCITY1, cr, constr);
987 wallcycle_stop(wcycle, ewcUPDATE);
988 constrain_velocities(step, nullptr, state, shake_vir, constr, bCalcVir, do_log, do_ene);
989 wallcycle_start(wcycle, ewcUPDATE);
990 /* if VV, compute the pressure and constraints */
991 /* For VV2, we strictly only need this if using pressure
992 * control, but we really would like to have accurate pressures
994 * Think about ways around this in the future?
995 * For now, keep this choice in comments.
997 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
998 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1000 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1001 if (bCalcEner && ir->eI == eiVVAK)
1003 bSumEkinhOld = TRUE;
1005 /* for vv, the first half of the integration actually corresponds to the previous step.
1006 So we need information from the last step in the first half of the integration */
1007 if (bGStat || do_per_step(step - 1, nstglobalcomm))
1009 wallcycle_stop(wcycle, ewcUPDATE);
1010 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1011 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1012 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller,
1013 state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
1014 (bGStat ? CGLO_GSTAT : 0) | (bCalcEner ? CGLO_ENERGY : 0)
1015 | (bTemp ? CGLO_TEMPERATURE : 0) | (bPres ? CGLO_PRESSURE : 0)
1016 | (bPres ? CGLO_CONSTRAINT : 0) | (bStopCM ? CGLO_STOPCM : 0)
1017 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1020 /* explanation of above:
1021 a) We compute Ekin at the full time step
1022 if 1) we are using the AveVel Ekin, and it's not the
1023 initial step, or 2) if we are using AveEkin, but need the full
1024 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1025 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1026 EkinAveVel because it's needed for the pressure */
1027 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1028 top_global, &top, state->x.rvec_array(), state->box,
1029 &shouldCheckNumberOfBondedInteractions);
1032 process_and_stopcm_grp(fplog, &vcm, *mdatoms, state->x.rvec_array(),
1033 state->v.rvec_array());
1034 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1036 wallcycle_start(wcycle, ewcUPDATE);
1038 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1043 m_add(force_vir, shake_vir,
1044 total_vir); /* we need the un-dispersion corrected total vir here */
1045 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ,
1046 trotter_seq, ettTSEQ2);
1048 /* TODO This is only needed when we're about to write
1049 * a checkpoint, because we use it after the restart
1050 * (in a kludge?). But what should we be doing if
1051 * the startingBehavior is NewSimulation or bInitStep are true? */
1052 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1054 copy_mat(shake_vir, state->svir_prev);
1055 copy_mat(force_vir, state->fvir_prev);
1057 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1059 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1060 enerd->term[F_TEMP] =
1061 sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1062 enerd->term[F_EKIN] = trace(ekind->ekin);
1065 else if (bExchanged)
1067 wallcycle_stop(wcycle, ewcUPDATE);
1068 /* We need the kinetic energy at minus the half step for determining
1069 * the full step kinetic energy and possibly for T-coupling.*/
1070 /* This may not be quite working correctly yet . . . . */
1071 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(),
1072 state->v.rvec_array(), state->box, state->lambda[efptVDW],
1073 mdatoms, nrnb, &vcm, wcycle, enerd, nullptr, nullptr, nullptr,
1074 nullptr, mu_tot, constr, &nullSignaller, state->box, nullptr,
1075 &bSumEkinhOld, CGLO_GSTAT | CGLO_TEMPERATURE);
1076 wallcycle_start(wcycle, ewcUPDATE);
1079 /* if it's the initial step, we performed this first step just to get the constraint virial */
1080 if (ir->eI == eiVV && bInitStep)
1082 copy_rvecn(vbuf, state->v.rvec_array(), 0, state->natoms);
1085 wallcycle_stop(wcycle, ewcUPDATE);
1088 /* compute the conserved quantity */
1091 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1094 last_ekin = enerd->term[F_EKIN];
1096 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1098 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1100 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1101 if (ir->efep != efepNO)
1103 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1107 /* ######## END FIRST UPDATE STEP ############## */
1108 /* ######## If doing VV, we now have v(dt) ###### */
1111 /* perform extended ensemble sampling in lambda - we don't
1112 actually move to the new state before outputting
1113 statistics, but if performing simulated tempering, we
1114 do update the velocities and the tau_t. */
1116 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state,
1117 state->dfhist, step, state->v.rvec_array(), mdatoms);
1118 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1121 copy_df_history(state_global->dfhist, state->dfhist);
1125 // Copy coordinate from the GPU for the output/checkpointing if the update is offloaded and
1126 // coordinates have not already been copied for i) search or ii) CPU force tasks.
1127 if (useGpuForUpdate && !bNS && !runScheduleWork->domainWork.haveCpuLocalForceWork
1128 && (do_per_step(step, ir->nstxout) || do_per_step(step, ir->nstxout_compressed)
1129 || checkpointHandler->isCheckpointingStep()))
1131 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1132 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1134 // Copy velocities if needed for the output/checkpointing.
1135 // NOTE: Copy on the search steps is done at the beginning of the step.
1136 if (useGpuForUpdate && !bNS
1137 && (do_per_step(step, ir->nstvout) || checkpointHandler->isCheckpointingStep()))
1139 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1140 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1142 // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
1143 // and update is offloaded hence forces are kept on the GPU for update and have not been
1144 // already transferred in do_force().
1145 // TODO: There should be an improved, explicit mechanism that ensures this copy is only executed
1146 // when the forces are ready on the GPU -- the same synchronizer should be used as the one
1147 // prior to GPU update.
1148 // TODO: When the output flags will be included in step workload, this copy can be combined with the
1149 // copy call in do_force(...).
1150 // NOTE: The forces should not be copied here if the vsites are present, since they were modified
1151 // on host after the D2H copy in do_force(...).
1152 if (runScheduleWork->stepWork.useGpuFBufferOps && (simulationWork.useGpuUpdate && !vsite)
1153 && do_per_step(step, ir->nstfout))
1155 stateGpu->copyForcesFromGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1156 stateGpu->waitForcesReadyOnHost(AtomLocality::Local);
1158 /* Now we have the energies and forces corresponding to the
1159 * coordinates at time t. We must output all of this before
1162 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state, state_global,
1163 observablesHistory, top_global, fr, outf, energyOutput, ekind, f,
1164 checkpointHandler->isCheckpointingStep(), bRerunMD, bLastStep,
1165 mdrunOptions.writeConfout, bSumEkinhOld);
1166 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1167 bInteractiveMDstep = imdSession->run(step, bNS, state->box, state->x.rvec_array(), t);
1169 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1170 if (startingBehavior != StartingBehavior::NewSimulation
1171 && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1173 copy_mat(state->svir_prev, shake_vir);
1174 copy_mat(state->fvir_prev, force_vir);
1177 stopHandler->setSignal();
1178 resetHandler->setSignal(walltime_accounting);
1180 if (bGStat || !PAR(cr))
1182 /* In parallel we only have to check for checkpointing in steps
1183 * where we do global communication,
1184 * otherwise the other nodes don't know.
1186 checkpointHandler->setSignal(walltime_accounting);
1189 /* ######### START SECOND UPDATE STEP ################# */
1191 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen
1192 controlled in preprocessing */
1194 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1196 gmx_bool bIfRandomize;
1197 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state->v, &upd, constr);
1198 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1199 if (constr && bIfRandomize)
1201 constrain_velocities(step, nullptr, state, tmp_vir, constr, bCalcVir, do_log, do_ene);
1204 /* Box is changed in update() when we do pressure coupling,
1205 * but we should still use the old box for energy corrections and when
1206 * writing it to the energy file, so it matches the trajectory files for
1207 * the same timestep above. Make a copy in a separate array.
1209 copy_mat(state->box, lastbox);
1213 wallcycle_start(wcycle, ewcUPDATE);
1214 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1217 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1218 /* We can only do Berendsen coupling after we have summed
1219 * the kinetic energy or virial. Since the happens
1220 * in global_state after update, we should only do it at
1221 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1226 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1227 update_pcouple_before_coordinates(fplog, step, ir, state, pressureCouplingMu, M, bInitStep);
1232 /* velocity half-step update */
1233 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1234 etrtVELOCITY2, cr, constr);
1237 /* Above, initialize just copies ekinh into ekin,
1238 * it doesn't copy position (for VV),
1239 * and entire integrator for MD.
1242 if (ir->eI == eiVVAK)
1244 cbuf.resize(state->x.size());
1245 std::copy(state->x.begin(), state->x.end(), cbuf.begin());
1248 /* With leap-frog type integrators we compute the kinetic energy
1249 * at a whole time step as the average of the half-time step kinetic
1250 * energies of two subsequent steps. Therefore we need to compute the
1251 * half step kinetic energy also if we need energies at the next step.
1253 const bool needHalfStepKineticEnergy =
1254 (!EI_VV(ir->eI) && (do_per_step(step + 1, nstglobalcomm) || step_rel + 1 == ir->nsteps));
1256 // Parrinello-Rahman requires the pressure to be availible before the update to compute
1257 // the velocity scaling matrix. Hence, it runs one step after the nstpcouple step.
1258 const bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN
1259 && do_per_step(step + ir->nstpcouple - 1, ir->nstpcouple));
1261 if (useGpuForUpdate)
1265 integrator->set(stateGpu->getCoordinates(), stateGpu->getVelocities(),
1266 stateGpu->getForces(), top.idef, *mdatoms, ekind->ngtc);
1268 // Copy data to the GPU after buffers might have being reinitialized
1269 stateGpu->copyVelocitiesToGpu(state->v, AtomLocality::Local);
1270 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1273 // If the buffer ops were not offloaded this step, the forces are on the host and have to be copied
1274 if (!runScheduleWork->stepWork.useGpuFBufferOps)
1276 stateGpu->copyForcesToGpu(ArrayRef<RVec>(f), AtomLocality::Local);
1279 const bool doTemperatureScaling =
1280 (ir->etc != etcNO && do_per_step(step + ir->nsttcouple - 1, ir->nsttcouple));
1282 // This applies Leap-Frog, LINCS and SETTLE in succession
1283 integrator->integrate(stateGpu->getForcesReadyOnDeviceEvent(
1284 AtomLocality::Local, runScheduleWork->stepWork.useGpuFBufferOps),
1285 ir->delta_t, true, bCalcVir, shake_vir, doTemperatureScaling,
1286 ekind->tcstat, doParrinelloRahman, ir->nstpcouple * ir->delta_t, M);
1288 // Copy velocities D2H after update if:
1289 // - Globals are computed this step (includes the energy output steps).
1290 // - Temperature is needed for the next step.
1291 if (bGStat || needHalfStepKineticEnergy)
1293 stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
1294 stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
1299 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1300 etrtPOSITION, cr, constr);
1302 wallcycle_stop(wcycle, ewcUPDATE);
1304 constrain_coordinates(step, &dvdl_constr, state, shake_vir, &upd, constr, bCalcVir,
1307 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state, cr, nrnb, wcycle, &upd,
1308 constr, do_log, do_ene);
1309 finish_update(ir, mdatoms, state, graph, nrnb, wcycle, &upd, constr);
1312 if (ir->bPull && ir->pull->bSetPbcRefToPrevStepCOM)
1314 updatePrevStepPullCom(pull_work, state);
1317 if (ir->eI == eiVVAK)
1319 /* erase F_EKIN and F_TEMP here? */
1320 /* just compute the kinetic energy at the half step to perform a trotter step */
1321 compute_globals(gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1322 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1323 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &nullSignaller, lastbox,
1324 nullptr, &bSumEkinhOld, (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE);
1325 wallcycle_start(wcycle, ewcUPDATE);
1326 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1327 /* now we know the scaling, we can compute the positions again */
1328 std::copy(cbuf.begin(), cbuf.end(), state->x.begin());
1330 update_coords(step, ir, mdatoms, state, f.arrayRefWithPadding(), fcd, ekind, M, &upd,
1331 etrtPOSITION, cr, constr);
1332 wallcycle_stop(wcycle, ewcUPDATE);
1334 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1335 /* are the small terms in the shake_vir here due
1336 * to numerical errors, or are they important
1337 * physically? I'm thinking they are just errors, but not completely sure.
1338 * For now, will call without actually constraining, constr=NULL*/
1339 finish_update(ir, mdatoms, state, graph, nrnb, wcycle, &upd, nullptr);
1343 /* this factor or 2 correction is necessary
1344 because half of the constraint force is removed
1345 in the vv step, so we have to double it. See
1346 the Redmine issue #1255. It is not yet clear
1347 if the factor of 2 is exact, or just a very
1348 good approximation, and this will be
1349 investigated. The next step is to see if this
1350 can be done adding a dhdl contribution from the
1351 rattle step, but this is somewhat more
1352 complicated with the current code. Will be
1353 investigated, hopefully for 4.6.3. However,
1354 this current solution is much better than
1355 having it completely wrong.
1357 enerd->term[F_DVDL_CONSTR] += 2 * dvdl_constr;
1361 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1364 if (vsite != nullptr)
1366 wallcycle_start(wcycle, ewcVSITECONSTR);
1367 if (graph != nullptr)
1369 shift_self(graph, state->box, state->x.rvec_array());
1371 construct_vsites(vsite, state->x.rvec_array(), ir->delta_t, state->v.rvec_array(),
1372 top.idef.iparams, top.idef.il, fr->ePBC, fr->bMolPBC, cr, state->box);
1374 if (graph != nullptr)
1376 unshift_self(graph, state->box, state->x.rvec_array());
1378 wallcycle_stop(wcycle, ewcVSITECONSTR);
1381 /* ############## IF NOT VV, Calculate globals HERE ############ */
1382 /* With Leap-Frog we can skip compute_globals at
1383 * non-communication steps, but we need to calculate
1384 * the kinetic energy one step before communication.
1387 // Organize to do inter-simulation signalling on steps if
1388 // and when algorithms require it.
1389 const bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1391 if (bGStat || needHalfStepKineticEnergy || doInterSimSignal)
1393 // Copy coordinates when needed to stop the CM motion.
1394 if (useGpuForUpdate && !EI_VV(ir->eI) && bStopCM)
1396 stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
1397 stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
1399 // Since we're already communicating at this step, we
1400 // can propagate intra-simulation signals. Note that
1401 // check_nstglobalcomm has the responsibility for
1402 // choosing the value of nstglobalcomm that is one way
1403 // bGStat becomes true, so we can't get into a
1404 // situation where e.g. checkpointing can't be
1406 bool doIntraSimSignal = true;
1407 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1410 gstat, cr, ir, fr, ekind, state->x.rvec_array(), state->v.rvec_array(),
1411 state->box, state->lambda[efptVDW], mdatoms, nrnb, &vcm, wcycle, enerd,
1412 force_vir, shake_vir, total_vir, pres, mu_tot, constr, &signaller, lastbox,
1413 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1414 (bGStat ? CGLO_GSTAT : 0) | (!EI_VV(ir->eI) && bCalcEner ? CGLO_ENERGY : 0)
1415 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1416 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1417 | (!EI_VV(ir->eI) ? CGLO_PRESSURE : 0) | CGLO_CONSTRAINT
1418 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
1420 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
1421 top_global, &top, state->x.rvec_array(), state->box,
1422 &shouldCheckNumberOfBondedInteractions);
1423 if (!EI_VV(ir->eI) && bStopCM)
1425 process_and_stopcm_grp(fplog, &vcm, *mdatoms, state->x.rvec_array(),
1426 state->v.rvec_array());
1427 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
1429 // TODO: The special case of removing CM motion should be dealt more gracefully
1430 if (useGpuForUpdate)
1432 stateGpu->copyCoordinatesToGpu(state->x, AtomLocality::Local);
1433 stateGpu->waitCoordinatesCopiedToDevice(AtomLocality::Local);
1439 /* ############# END CALC EKIN AND PRESSURE ################# */
1441 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1442 the virial that should probably be addressed eventually. state->veta has better properies,
1443 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1444 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1446 if (ir->efep != efepNO && !EI_VV(ir->eI))
1448 /* Sum up the foreign energy and dhdl terms for md and sd.
1449 Currently done every step so that dhdl is correct in the .edr */
1450 sum_dhdl(enerd, state->lambda, *ir->fepvals);
1453 update_pcouple_after_coordinates(fplog, step, ir, mdatoms, pres, force_vir, shake_vir,
1454 pressureCouplingMu, state, nrnb, &upd, !useGpuForUpdate);
1456 const bool doBerendsenPressureCoupling =
1457 (inputrec->epc == epcBERENDSEN && do_per_step(step, inputrec->nstpcouple));
1458 if (useGpuForUpdate && (doBerendsenPressureCoupling || doParrinelloRahman))
1460 integrator->scaleCoordinates(pressureCouplingMu);
1462 set_pbc(&pbc, epbcXYZ, state->box);
1463 integrator->setPbc(&pbc);
1466 /* ################# END UPDATE STEP 2 ################# */
1467 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1469 /* The coordinates (x) were unshifted in update */
1472 /* We will not sum ekinh_old,
1473 * so signal that we still have to do it.
1475 bSumEkinhOld = TRUE;
1480 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1482 /* use the directly determined last velocity, not actually the averaged half steps */
1483 if (bTrotter && ir->eI == eiVV)
1485 enerd->term[F_EKIN] = last_ekin;
1487 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1489 if (integratorHasConservedEnergyQuantity(ir))
1493 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1497 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1500 /* ######### END PREPARING EDR OUTPUT ########### */
1506 if (fplog && do_log && bDoExpanded)
1508 /* only needed if doing expanded ensemble */
1509 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals,
1510 ir->bSimTemp ? ir->simtempvals : nullptr,
1511 state_global->dfhist, state->fep_state, ir->nstlog, step);
1515 energyOutput.addDataAtEnergyStep(bDoDHDL, bCalcEnerStep, t, mdatoms->tmass, enerd, state,
1516 ir->fepvals, ir->expandedvals, lastbox, shake_vir,
1517 force_vir, total_vir, pres, ekind, mu_tot, constr);
1521 energyOutput.recordNonEnergyStep();
1524 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1525 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1527 if (doSimulatedAnnealing)
1529 energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
1531 if (do_log || do_ene || do_dr || do_or)
1533 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
1534 do_log ? fplog : nullptr, step, t, fcd, awh.get());
1539 pull_print_output(pull_work, step, t);
1542 if (do_per_step(step, ir->nstlog))
1544 if (fflush(fplog) != 0)
1546 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1552 /* Have to do this part _after_ outputting the logfile and the edr file */
1553 /* Gets written into the state at the beginning of next loop*/
1554 state->fep_state = lamnew;
1556 /* Print the remaining wall clock time for the run */
1557 if (isMasterSimMasterRank(ms, MASTER(cr)) && (do_verbose || gmx_got_usr_signal()) && !bPMETunePrinting)
1561 fprintf(stderr, "\n");
1563 print_time(stderr, walltime_accounting, step, ir, cr);
1566 /* Ion/water position swapping.
1567 * Not done in last step since trajectory writing happens before this call
1568 * in the MD loop and exchanges would be lost anyway. */
1569 bNeedRepartition = FALSE;
1570 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep && do_per_step(step, ir->swap->nstswap))
1573 do_swapcoords(cr, step, t, ir, swap, wcycle, as_rvec_array(state->x.data()),
1574 state->box, MASTER(cr) && mdrunOptions.verbose, bRerunMD);
1576 if (bNeedRepartition && DOMAINDECOMP(cr))
1578 dd_collect_state(cr->dd, state, state_global);
1582 /* Replica exchange */
1586 bExchanged = replica_exchange(fplog, cr, ms, repl_ex, state_global, enerd, state, step, t);
1589 if ((bExchanged || bNeedRepartition) && DOMAINDECOMP(cr))
1591 dd_partition_system(fplog, mdlog, step, cr, TRUE, 1, state_global, *top_global, ir,
1592 imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
1593 nrnb, wcycle, FALSE);
1594 shouldCheckNumberOfBondedInteractions = true;
1595 upd.setNumAtoms(state->natoms);
1601 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1602 /* With all integrators, except VV, we need to retain the pressure
1603 * at the current step for coupling at the next step.
1605 if ((state->flags & (1U << estPRES_PREV))
1606 && (bGStatEveryStep || (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1608 /* Store the pressure in t_state for pressure coupling
1609 * at the next MD step.
1611 copy_mat(pres, state->pres_prev);
1614 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1616 if ((membed != nullptr) && (!bLastStep))
1618 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1621 cycles = wallcycle_stop(wcycle, ewcSTEP);
1622 if (DOMAINDECOMP(cr) && wcycle)
1624 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1627 /* increase the MD step number */
1631 resetHandler->resetCounters(step, step_rel, mdlog, fplog, cr, fr->nbv.get(), nrnb,
1632 fr->pmedata, pme_loadbal, wcycle, walltime_accounting);
1634 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1635 imdSession->updateEnergyRecordAndSendPositionsAndEnergies(bInteractiveMDstep, step, bCalcEner);
1637 /* End of main MD loop */
1639 /* Closing TNG files can include compressing data. Therefore it is good to do that
1640 * before stopping the time measurements. */
1641 mdoutf_tng_close(outf);
1643 /* Stop measuring walltime */
1644 walltime_accounting_end_time(walltime_accounting);
1646 if (!thisRankHasDuty(cr, DUTY_PME))
1648 /* Tell the PME only node to finish */
1649 gmx_pme_send_finish(cr);
1654 if (ir->nstcalcenergy > 0)
1656 energyOutput.printAnnealingTemperatures(fplog, groups, &(ir->opts));
1657 energyOutput.printAverages(fplog, groups);
1664 pme_loadbal_done(pme_loadbal, fplog, mdlog, fr->nbv->useGpu());
1667 done_shellfc(fplog, shellfc, step_rel);
1669 if (useReplicaExchange && MASTER(cr))
1671 print_replica_exchange_statistics(fplog, repl_ex);
1674 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1676 global_stat_destroy(gstat);