2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * \brief Implements the loop for simulation reruns
39 * \author Pascal Merz <pascal.merz@colorado.edu>
40 * \ingroup module_mdrun
52 #include "gromacs/applied_forces/awh/awh.h"
53 #include "gromacs/commandline/filenm.h"
54 #include "gromacs/domdec/collect.h"
55 #include "gromacs/domdec/dlbtiming.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/domdec/mdsetup.h"
60 #include "gromacs/domdec/partition.h"
61 #include "gromacs/essentialdynamics/edsam.h"
62 #include "gromacs/ewald/pme_load_balancing.h"
63 #include "gromacs/ewald/pme_pp.h"
64 #include "gromacs/fileio/trxio.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gmxlib/nrnb.h"
67 #include "gromacs/gpu_utils/gpu_utils.h"
68 #include "gromacs/listed_forces/listed_forces.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/utilities.h"
71 #include "gromacs/math/vec.h"
72 #include "gromacs/math/vectypes.h"
73 #include "gromacs/mdlib/checkpointhandler.h"
74 #include "gromacs/mdlib/compute_io.h"
75 #include "gromacs/mdlib/constr.h"
76 #include "gromacs/mdlib/ebin.h"
77 #include "gromacs/mdlib/enerdata_utils.h"
78 #include "gromacs/mdlib/energyoutput.h"
79 #include "gromacs/mdlib/expanded.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/freeenergyparameters.h"
84 #include "gromacs/mdlib/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/mdoutf.h"
87 #include "gromacs/mdlib/membed.h"
88 #include "gromacs/mdlib/resethandler.h"
89 #include "gromacs/mdlib/sighandler.h"
90 #include "gromacs/mdlib/simulationsignal.h"
91 #include "gromacs/mdlib/stat.h"
92 #include "gromacs/mdlib/stophandler.h"
93 #include "gromacs/mdlib/tgroup.h"
94 #include "gromacs/mdlib/trajectory_writing.h"
95 #include "gromacs/mdlib/update.h"
96 #include "gromacs/mdlib/vcm.h"
97 #include "gromacs/mdlib/vsite.h"
98 #include "gromacs/mdrunutility/handlerestart.h"
99 #include "gromacs/mdrunutility/multisim.h"
100 #include "gromacs/mdrunutility/printtime.h"
101 #include "gromacs/mdtypes/awh_history.h"
102 #include "gromacs/mdtypes/awh_params.h"
103 #include "gromacs/mdtypes/commrec.h"
104 #include "gromacs/mdtypes/df_history.h"
105 #include "gromacs/mdtypes/energyhistory.h"
106 #include "gromacs/mdtypes/forcebuffers.h"
107 #include "gromacs/mdtypes/forcerec.h"
108 #include "gromacs/mdtypes/group.h"
109 #include "gromacs/mdtypes/inputrec.h"
110 #include "gromacs/mdtypes/interaction_const.h"
111 #include "gromacs/mdtypes/md_enums.h"
112 #include "gromacs/mdtypes/mdatom.h"
113 #include "gromacs/mdtypes/mdrunoptions.h"
114 #include "gromacs/mdtypes/observableshistory.h"
115 #include "gromacs/mdtypes/simulation_workload.h"
116 #include "gromacs/mdtypes/state.h"
117 #include "gromacs/mimic/utilities.h"
118 #include "gromacs/pbcutil/pbc.h"
119 #include "gromacs/pulling/pull.h"
120 #include "gromacs/swap/swapcoords.h"
121 #include "gromacs/timing/wallcycle.h"
122 #include "gromacs/timing/walltime_accounting.h"
123 #include "gromacs/topology/atoms.h"
124 #include "gromacs/topology/idef.h"
125 #include "gromacs/topology/mtop_util.h"
126 #include "gromacs/topology/topology.h"
127 #include "gromacs/trajectory/trajectoryframe.h"
128 #include "gromacs/utility/basedefinitions.h"
129 #include "gromacs/utility/cstringutil.h"
130 #include "gromacs/utility/fatalerror.h"
131 #include "gromacs/utility/logger.h"
132 #include "gromacs/utility/real.h"
134 #include "legacysimulator.h"
135 #include "replicaexchange.h"
138 using gmx::SimulationSignaller;
139 using gmx::VirtualSitesHandler;
141 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
143 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
144 * \param[in,out] globalState The global state container
145 * \param[in] constructVsites When true, vsite coordinates are constructed
146 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
148 static void prepareRerunState(const t_trxframe& rerunFrame,
149 t_state* globalState,
150 bool constructVsites,
151 const VirtualSitesHandler* vsite)
153 auto x = makeArrayRef(globalState->x);
154 auto rerunX = arrayRefFromArray(reinterpret_cast<gmx::RVec*>(rerunFrame.x), globalState->natoms);
155 std::copy(rerunX.begin(), rerunX.end(), x.begin());
156 copy_mat(rerunFrame.box, globalState->box);
160 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
162 vsite->construct(globalState->x, globalState->v, globalState->box, gmx::VSiteOperation::PositionsAndVelocities);
166 void gmx::LegacySimulator::do_rerun()
168 // TODO Historically, the EM and MD "integrators" used different
169 // names for the t_inputrec *parameter, but these must have the
170 // same name, now that it's a member of a struct. We use this ir
171 // alias to avoid a large ripple of nearly useless changes.
172 // t_inputrec is being replaced by IMdpOptionsProvider, so this
173 // will go away eventually.
174 const t_inputrec* ir = inputrec;
175 int64_t step, step_rel;
177 bool isLastStep = false;
178 bool doFreeEnergyPerturbation = false;
179 unsigned int force_flags;
180 tensor force_vir, shake_vir, total_vir, pres;
181 t_trxstatus* status = nullptr;
184 gmx_localtop_t top(top_global.ffparams);
186 gmx_global_stat_t gstat;
187 gmx_shellfc_t* shellfc;
191 SimulationSignals signals;
192 // Most global communnication stages don't propagate mdrun
193 // signals, and will use this object to achieve that.
194 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
199 "Note that it is planned that the command gmx mdrun -rerun will "
200 "be available in a different form in a future version of GROMACS, "
201 "e.g. gmx rerun -f.");
203 if (ir->efep != FreeEnergyPerturbationType::No
204 && (mdAtoms->mdatoms()->nMassPerturbed > 0 || (constr && constr->havePerturbedConstraints())))
207 "Perturbed masses or constraints are not supported by rerun. "
208 "Either make a .tpr without mass and constraint perturbation, "
209 "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
213 gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
217 gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
221 gmx_fatal(FARGS, "AWH not supported by rerun.");
223 if (replExParams.exchangeInterval > 0)
225 gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
227 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
229 gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
233 gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
237 gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
239 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc, [](SimulatedAnnealing i) {
240 return i != SimulatedAnnealing::No;
243 gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
246 /* Rerun can't work if an output file name is the same as the input file name.
247 * If this is the case, the user will get an error telling them what the issue is.
249 if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0
250 || strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
253 "When using mdrun -rerun, the name of the input trajectory file "
254 "%s cannot be identical to the name of an output file (whether "
255 "given explicitly with -o or -x, or by default)",
256 opt2fn("-rerun", nfile, fnm));
259 /* Settings for rerun */
261 // TODO: Avoid changing inputrec (#3854)
262 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
263 nonConstInputrec->nstlist = 1;
264 nonConstInputrec->nstcalcenergy = 1;
265 nonConstInputrec->nstxout_compressed = 0;
267 int nstglobalcomm = 1;
268 const bool bNS = true;
270 const SimulationGroups* groups = &top_global.groups;
271 if (ir->eI == IntegrationAlgorithm::Mimic)
273 auto nonConstGlobalTopology = const_cast<gmx_mtop_t*>(&top_global);
274 nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(top_global);
276 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
277 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
278 initialize_lambdas(fplog, *ir, MASTER(cr), fep_state, lambda);
279 const bool simulationsShareState = false;
280 gmx_mdoutf* outf = init_mdoutf(fplog,
291 StartingBehavior::NewSimulation,
292 simulationsShareState,
294 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
298 mdoutf_get_fp_dhdl(outf),
300 StartingBehavior::NewSimulation,
301 simulationsShareState,
304 gstat = global_stat_init(ir);
306 /* Check for polarizable models and flexible constraints */
307 shellfc = init_shell_flexcon(fplog,
309 constr ? constr->numFlexibleConstraints() : 0,
312 runScheduleWork->simulationWork.useGpuPme);
315 double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
316 if ((io > 2000) && MASTER(cr))
318 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
322 // Local state only becomes valid now.
323 std::unique_ptr<t_state> stateInstance;
326 if (DOMAINDECOMP(cr))
328 stateInstance = std::make_unique<t_state>();
329 state = stateInstance.get();
330 dd_init_local_state(*cr->dd, state_global, state);
332 /* Distribute the charge groups over the nodes from the master node */
333 dd_partition_system(fplog,
357 state_change_natoms(state_global, state_global->natoms);
358 /* Copy the pointer to the global state */
359 state = state_global;
361 mdAlgorithmsSetupAtomData(cr, *ir, top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
364 auto mdatoms = mdAtoms->mdatoms();
366 // NOTE: The global state is no longer used at this point.
367 // But state_global is still used as temporary storage space for writing
368 // the global state to file and potentially for replica exchange.
369 // (Global topology should persist.)
371 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
373 if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
375 doFreeEnergyPerturbation = true;
379 int cglo_flags = CGLO_GSTAT;
380 if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
382 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
384 bool bSumEkinhOld = false;
385 t_vcm* vcm = nullptr;
386 compute_globals(gstat,
391 makeConstArrayRef(state->x),
392 makeConstArrayRef(state->v),
403 gmx::ArrayRef<real>{},
408 if (DOMAINDECOMP(cr))
410 checkNumberOfBondedInteractions(
411 mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
418 "starting md rerun '%s', reading coordinates from"
419 " input trajectory '%s'\n\n",
421 opt2fn("-rerun", nfile, fnm));
422 if (mdrunOptions.verbose)
425 "Calculated time to finish depends on nsteps from "
426 "run input file,\nwhich may not correspond to the time "
427 "needed to process input trajectory.\n\n");
429 fprintf(fplog, "\n");
432 walltime_accounting_start_time(walltime_accounting);
433 wallcycle_start(wcycle, ewcRUN);
434 print_start(fplog, cr, walltime_accounting, "mdrun");
436 /***********************************************************
440 ************************************************************/
446 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
452 isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
453 if (rerun_fr.natoms != top_global.natoms)
456 "Number of atoms in trajectory (%d) does not match the "
457 "run input file (%d)\n",
462 if (ir->pbcType != PbcType::No)
467 "Rerun trajectory frame step %" PRId64
469 "does not contain a box, while pbc is used",
473 if (max_cutoff2(ir->pbcType, rerun_fr.box) < gmx::square(fr->rlist))
476 "Rerun trajectory frame step %" PRId64
478 "has too small box dimensions",
488 "Rerun does not report kinetic energy, total energy, temperature, virial and "
493 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
496 if (ir->pbcType != PbcType::No)
498 /* Set the shift vectors.
499 * Necessary here when have a static box different from the tpr box.
501 calc_shifts(rerun_fr.box, fr->shift_vec);
504 step = ir->init_step;
507 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
508 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
512 mdrunOptions.reproducible,
514 mdrunOptions.maximumHoursToRun,
519 walltime_accounting);
521 // we don't do counter resetting in rerun - finish will always be valid
522 walltime_accounting_set_valid_finish(walltime_accounting);
524 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
526 /* and stop now if we should */
527 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
530 wallcycle_start(wcycle, ewcSTEP);
534 step = rerun_fr.step;
535 step_rel = step - ir->init_step;
546 if (ir->efep != FreeEnergyPerturbationType::No && MASTER(cr))
548 if (rerun_fr.bLambda)
550 ir->fepvals->init_lambda = rerun_fr.lambda;
554 if (rerun_fr.bFepState)
556 state->fep_state = rerun_fr.fep_state;
560 state_global->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
565 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
566 if (constructVsites && DOMAINDECOMP(cr))
569 "Vsite recalculation with -rerun is not implemented with domain "
571 "use a single rank");
573 prepareRerunState(rerun_fr, state_global, constructVsites, vsite);
576 isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
578 if (DOMAINDECOMP(cr))
580 /* Repartition the domain decomposition */
581 const bool bMasterState = true;
582 dd_partition_system(fplog,
602 mdrunOptions.verbose);
607 EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
610 if (ir->efep != FreeEnergyPerturbationType::No)
612 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
615 force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
616 | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
617 GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
621 /* Now is the time to relax the shells */
622 relax_shell_flexcon(fplog,
625 mdrunOptions.verbose,
637 state->x.arrayRefWithPadding(),
638 state->v.arrayRefWithPadding(),
653 ddBalanceRegionHandler);
657 /* The coordinates (x) are shifted (to get whole molecules)
659 * This is parallellized as well, and does communication too.
660 * Check comments in sim_util.c
663 gmx_edsam* ed = nullptr;
677 state->x.arrayRefWithPadding(),
690 GMX_FORCE_NS | force_flags,
691 ddBalanceRegionHandler);
694 /* Now we have the energies and forces corresponding to the
695 * coordinates at time t.
698 const bool isCheckpointingStep = false;
699 const bool doRerun = true;
700 const bool bSumEkinhOld = false;
701 do_md_trajectory_writing(fplog,
721 mdrunOptions.writeConfout,
725 stopHandler->setSignal();
728 const bool doInterSimSignal = false;
729 const bool doIntraSimSignal = true;
730 bool bSumEkinhOld = false;
731 t_vcm* vcm = nullptr;
732 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
734 int cglo_flags = CGLO_GSTAT | CGLO_ENERGY;
735 if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
737 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
739 compute_globals(gstat,
744 makeConstArrayRef(state->x),
745 makeConstArrayRef(state->v),
756 constr != nullptr ? constr->rmsdData() : gmx::ArrayRef<real>{},
761 if (DOMAINDECOMP(cr))
763 checkNumberOfBondedInteractions(
764 mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
768 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
769 the virial that should probably be addressed eventually. state->veta has better properies,
770 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
771 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
776 const bool bCalcEnerStep = true;
777 energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
783 ir->expandedvals.get(),
785 PTCouplingArrays({ state->boxv,
786 state->nosehoover_xi,
787 state->nosehoover_vxi,
789 state->nhpres_vxi }),
799 const bool do_ene = true;
800 const bool do_log = true;
802 const bool do_dr = ir->nstdisreout != 0;
803 const bool do_or = ir->nstorireout != 0;
805 EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
806 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
810 do_log ? fplog : nullptr,
816 if (do_per_step(step, ir->nstlog))
818 if (fflush(fplog) != 0)
820 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
825 /* Print the remaining wall clock time for the run */
826 if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
830 fprintf(stderr, "\n");
832 print_time(stderr, walltime_accounting, step, ir, cr);
835 /* Ion/water position swapping.
836 * Not done in last step since trajectory writing happens before this call
837 * in the MD loop and exchanges would be lost anyway. */
838 if ((ir->eSwapCoords != SwapType::No) && (step > 0) && !isLastStep
839 && do_per_step(step, ir->swap->nstswap))
841 const bool doRerun = true;
850 MASTER(cr) && mdrunOptions.verbose,
856 /* read next frame from input trajectory */
857 isLastStep = !read_next_frame(oenv, status, &rerun_fr);
862 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
865 cycles = wallcycle_stop(wcycle, ewcSTEP);
866 if (DOMAINDECOMP(cr) && wcycle)
868 dd_cycles_add(cr->dd, cycles, ddCyclStep);
873 /* increase the MD step number */
878 /* End of main MD loop */
880 /* Closing TNG files can include compressing data. Therefore it is good to do that
881 * before stopping the time measurements. */
882 mdoutf_tng_close(outf);
884 /* Stop measuring walltime */
885 walltime_accounting_end_time(walltime_accounting);
892 if (!thisRankHasDuty(cr, DUTY_PME))
894 /* Tell the PME only node to finish */
895 gmx_pme_send_finish(cr);
900 done_shellfc(fplog, shellfc, step_rel);
902 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);