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/localtopologychecker.h"
60 #include "gromacs/domdec/mdsetup.h"
61 #include "gromacs/domdec/partition.h"
62 #include "gromacs/essentialdynamics/edsam.h"
63 #include "gromacs/ewald/pme_load_balancing.h"
64 #include "gromacs/ewald/pme_pp.h"
65 #include "gromacs/fileio/trxio.h"
66 #include "gromacs/gmxlib/network.h"
67 #include "gromacs/gmxlib/nrnb.h"
68 #include "gromacs/gpu_utils/gpu_utils.h"
69 #include "gromacs/listed_forces/listed_forces.h"
70 #include "gromacs/math/functions.h"
71 #include "gromacs/math/utilities.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/math/vectypes.h"
74 #include "gromacs/mdlib/checkpointhandler.h"
75 #include "gromacs/mdlib/compute_io.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/ebin.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/energyoutput.h"
80 #include "gromacs/mdlib/expanded.h"
81 #include "gromacs/mdlib/force.h"
82 #include "gromacs/mdlib/force_flags.h"
83 #include "gromacs/mdlib/forcerec.h"
84 #include "gromacs/mdlib/freeenergyparameters.h"
85 #include "gromacs/mdlib/md_support.h"
86 #include "gromacs/mdlib/mdatoms.h"
87 #include "gromacs/mdlib/mdoutf.h"
88 #include "gromacs/mdlib/membed.h"
89 #include "gromacs/mdlib/resethandler.h"
90 #include "gromacs/mdlib/sighandler.h"
91 #include "gromacs/mdlib/simulationsignal.h"
92 #include "gromacs/mdlib/stat.h"
93 #include "gromacs/mdlib/stophandler.h"
94 #include "gromacs/mdlib/tgroup.h"
95 #include "gromacs/mdlib/trajectory_writing.h"
96 #include "gromacs/mdlib/update.h"
97 #include "gromacs/mdlib/vcm.h"
98 #include "gromacs/mdlib/vsite.h"
99 #include "gromacs/mdrunutility/handlerestart.h"
100 #include "gromacs/mdrunutility/multisim.h"
101 #include "gromacs/mdrunutility/printtime.h"
102 #include "gromacs/mdtypes/awh_history.h"
103 #include "gromacs/mdtypes/awh_params.h"
104 #include "gromacs/mdtypes/commrec.h"
105 #include "gromacs/mdtypes/df_history.h"
106 #include "gromacs/mdtypes/energyhistory.h"
107 #include "gromacs/mdtypes/forcebuffers.h"
108 #include "gromacs/mdtypes/forcerec.h"
109 #include "gromacs/mdtypes/group.h"
110 #include "gromacs/mdtypes/inputrec.h"
111 #include "gromacs/mdtypes/interaction_const.h"
112 #include "gromacs/mdtypes/md_enums.h"
113 #include "gromacs/mdtypes/mdatom.h"
114 #include "gromacs/mdtypes/mdrunoptions.h"
115 #include "gromacs/mdtypes/observableshistory.h"
116 #include "gromacs/mdtypes/observablesreducer.h"
117 #include "gromacs/mdtypes/simulation_workload.h"
118 #include "gromacs/mdtypes/state.h"
119 #include "gromacs/mimic/utilities.h"
120 #include "gromacs/pbcutil/pbc.h"
121 #include "gromacs/pulling/output.h"
122 #include "gromacs/pulling/pull.h"
123 #include "gromacs/swap/swapcoords.h"
124 #include "gromacs/timing/wallcycle.h"
125 #include "gromacs/timing/walltime_accounting.h"
126 #include "gromacs/topology/atoms.h"
127 #include "gromacs/topology/idef.h"
128 #include "gromacs/topology/mtop_util.h"
129 #include "gromacs/topology/topology.h"
130 #include "gromacs/trajectory/trajectoryframe.h"
131 #include "gromacs/utility/basedefinitions.h"
132 #include "gromacs/utility/cstringutil.h"
133 #include "gromacs/utility/fatalerror.h"
134 #include "gromacs/utility/logger.h"
135 #include "gromacs/utility/real.h"
137 #include "legacysimulator.h"
138 #include "replicaexchange.h"
141 using gmx::SimulationSignaller;
142 using gmx::VirtualSitesHandler;
144 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
146 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
147 * \param[in,out] globalState The global state container
148 * \param[in] constructVsites When true, vsite coordinates are constructed
149 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
151 static void prepareRerunState(const t_trxframe& rerunFrame,
152 t_state* globalState,
153 bool constructVsites,
154 const VirtualSitesHandler* vsite)
156 auto x = makeArrayRef(globalState->x);
157 auto rerunX = arrayRefFromArray(reinterpret_cast<gmx::RVec*>(rerunFrame.x), globalState->natoms);
158 std::copy(rerunX.begin(), rerunX.end(), x.begin());
159 copy_mat(rerunFrame.box, globalState->box);
163 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
165 vsite->construct(globalState->x, globalState->v, globalState->box, gmx::VSiteOperation::PositionsAndVelocities);
169 void gmx::LegacySimulator::do_rerun()
171 // TODO Historically, the EM and MD "integrators" used different
172 // names for the t_inputrec *parameter, but these must have the
173 // same name, now that it's a member of a struct. We use this ir
174 // alias to avoid a large ripple of nearly useless changes.
175 // t_inputrec is being replaced by IMdpOptionsProvider, so this
176 // will go away eventually.
177 const t_inputrec* ir = inputrec;
179 bool isLastStep = false;
180 bool doFreeEnergyPerturbation = false;
181 unsigned int force_flags;
182 tensor force_vir, shake_vir, total_vir, pres;
183 t_trxstatus* status = nullptr;
187 gmx_global_stat_t gstat;
188 gmx_shellfc_t* shellfc;
192 SimulationSignals signals;
193 // Most global communnication stages don't propagate mdrun
194 // signals, and will use this object to achieve that.
195 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
200 "Note that it is planned that the command gmx mdrun -rerun will "
201 "be available in a different form in a future version of GROMACS, "
202 "e.g. gmx rerun -f.");
204 if (ir->efep != FreeEnergyPerturbationType::No
205 && (mdAtoms->mdatoms()->nMassPerturbed > 0 || (constr && constr->havePerturbedConstraints())))
208 "Perturbed masses or constraints are not supported by rerun. "
209 "Either make a .tpr without mass and constraint perturbation, "
210 "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
214 gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
218 gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
222 gmx_fatal(FARGS, "AWH not supported by rerun.");
224 if (replExParams.exchangeInterval > 0)
226 gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
228 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
230 gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
234 gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
238 gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
240 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc, [](SimulatedAnnealing i) {
241 return i != SimulatedAnnealing::No;
244 gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
247 /* Rerun can't work if an output file name is the same as the input file name.
248 * If this is the case, the user will get an error telling them what the issue is.
250 if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0
251 || strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
254 "When using mdrun -rerun, the name of the input trajectory file "
255 "%s cannot be identical to the name of an output file (whether "
256 "given explicitly with -o or -x, or by default)",
257 opt2fn("-rerun", nfile, fnm));
260 /* Settings for rerun */
262 // TODO: Avoid changing inputrec (#3854)
263 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
264 nonConstInputrec->nstlist = 1;
265 nonConstInputrec->nstcalcenergy = 1;
266 nonConstInputrec->nstxout_compressed = 0;
268 int nstglobalcomm = 1;
269 const bool bNS = true;
271 ObservablesReducer observablesReducer = observablesReducerBuilder->build();
273 const SimulationGroups* groups = &top_global.groups;
274 if (ir->eI == IntegrationAlgorithm::Mimic)
276 auto* nonConstGlobalTopology = const_cast<gmx_mtop_t*>(&top_global);
277 nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(top_global);
279 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
280 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
281 initialize_lambdas(fplog,
285 ir->simtempvals->temperatures,
286 gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
290 const bool simulationsShareState = false;
291 gmx_mdoutf* outf = init_mdoutf(fplog,
302 StartingBehavior::NewSimulation,
303 simulationsShareState,
305 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
309 mdoutf_get_fp_dhdl(outf),
311 StartingBehavior::NewSimulation,
312 simulationsShareState,
315 gstat = global_stat_init(ir);
317 /* Check for polarizable models and flexible constraints */
318 shellfc = init_shell_flexcon(fplog,
320 constr ? constr->numFlexibleConstraints() : 0,
322 haveDDAtomOrdering(*cr),
323 runScheduleWork->simulationWork.useGpuPme);
326 double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
327 if ((io > 2000) && MASTER(cr))
329 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
333 if (haveDDAtomOrdering(*cr))
335 // Local state only becomes valid now.
336 dd_init_local_state(*cr->dd, state_global, state);
338 /* Distribute the charge groups over the nodes from the master node */
339 dd_partition_system(fplog,
363 state_change_natoms(state_global, state_global->natoms);
364 /* Copy the pointer to the global state */
365 state = state_global;
367 mdAlgorithmsSetupAtomData(cr, *ir, top_global, top, fr, &f, mdAtoms, constr, vsite, shellfc);
370 auto* mdatoms = mdAtoms->mdatoms();
371 fr->longRangeNonbondeds->updateAfterPartition(*mdatoms);
373 // NOTE: The global state is no longer used at this point.
374 // But state_global is still used as temporary storage space for writing
375 // the global state to file and potentially for replica exchange.
376 // (Global topology should persist.)
378 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
380 if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
382 doFreeEnergyPerturbation = true;
385 int64_t step = ir->init_step;
386 int64_t step_rel = 0;
389 int cglo_flags = CGLO_GSTAT;
390 bool bSumEkinhOld = false;
391 t_vcm* vcm = nullptr;
392 compute_globals(gstat,
397 makeConstArrayRef(state->x),
398 makeConstArrayRef(state->v),
414 &observablesReducer);
415 // Clean up after pre-step use of compute_globals()
416 observablesReducer.markAsReadyToReduce();
422 "starting md rerun '%s', reading coordinates from"
423 " input trajectory '%s'\n\n",
425 opt2fn("-rerun", nfile, fnm));
426 if (mdrunOptions.verbose)
429 "Calculated time to finish depends on nsteps from "
430 "run input file,\nwhich may not correspond to the time "
431 "needed to process input trajectory.\n\n");
433 fprintf(fplog, "\n");
436 walltime_accounting_start_time(walltime_accounting);
437 wallcycle_start(wcycle, WallCycleCounter::Run);
438 print_start(fplog, cr, walltime_accounting, "mdrun");
440 /***********************************************************
444 ************************************************************/
450 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
456 isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
457 if (rerun_fr.natoms != top_global.natoms)
460 "Number of atoms in trajectory (%d) does not match the "
461 "run input file (%d)\n",
466 if (ir->pbcType != PbcType::No)
471 "Rerun trajectory frame step %" PRId64
473 "does not contain a box, while pbc is used",
477 if (max_cutoff2(ir->pbcType, rerun_fr.box) < gmx::square(fr->rlist))
480 "Rerun trajectory frame step %" PRId64
482 "has too small box dimensions",
492 "Rerun does not report kinetic energy, total energy, temperature, virial and "
497 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
500 if (ir->pbcType != PbcType::No)
502 /* Set the shift vectors.
503 * Necessary here when have a static box different from the tpr box.
505 calc_shifts(rerun_fr.box, fr->shift_vec);
508 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
509 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
513 mdrunOptions.reproducible,
515 mdrunOptions.maximumHoursToRun,
520 walltime_accounting);
522 // we don't do counter resetting in rerun - finish will always be valid
523 walltime_accounting_set_valid_finish(walltime_accounting);
525 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
527 /* and stop now if we should */
528 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
531 wallcycle_start(wcycle, WallCycleCounter::Step);
535 step = rerun_fr.step;
536 step_rel = step - ir->init_step;
547 if (ir->efep != FreeEnergyPerturbationType::No && MASTER(cr))
549 if (rerun_fr.bLambda)
551 ir->fepvals->init_lambda = rerun_fr.lambda;
555 if (rerun_fr.bFepState)
557 state->fep_state = rerun_fr.fep_state;
561 state_global->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
566 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
567 if (constructVsites && haveDDAtomOrdering(*cr))
570 "Vsite recalculation with -rerun is not implemented with domain "
572 "use a single rank");
574 prepareRerunState(rerun_fr, state_global, constructVsites, vsite);
577 isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
579 if (haveDDAtomOrdering(*cr))
581 /* Repartition the domain decomposition */
582 const bool bMasterState = true;
583 dd_partition_system(fplog,
603 mdrunOptions.verbose);
608 EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
611 if (ir->efep != FreeEnergyPerturbationType::No)
613 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
616 fr->longRangeNonbondeds->updateAfterPartition(*mdatoms);
618 force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
619 | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
620 GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
624 /* Now is the time to relax the shells */
625 relax_shell_flexcon(fplog,
628 mdrunOptions.verbose,
640 state->x.arrayRefWithPadding(),
641 state->v.arrayRefWithPadding(),
648 fr->longRangeNonbondeds.get(),
657 ddBalanceRegionHandler);
661 /* The coordinates (x) are shifted (to get whole molecules)
663 * This is parallellized as well, and does communication too.
664 * Check comments in sim_util.c
667 gmx_edsam* ed = nullptr;
681 state->x.arrayRefWithPadding(),
694 fr->longRangeNonbondeds.get(),
695 GMX_FORCE_NS | force_flags,
696 ddBalanceRegionHandler);
699 /* Now we have the energies and forces corresponding to the
700 * coordinates at time t.
703 const bool isCheckpointingStep = false;
704 const bool doRerun = true;
705 const bool bSumEkinhOld = false;
706 do_md_trajectory_writing(fplog,
726 mdrunOptions.writeConfout,
730 stopHandler->setSignal();
733 const bool doInterSimSignal = false;
734 const bool doIntraSimSignal = true;
735 bool bSumEkinhOld = false;
736 t_vcm* vcm = nullptr;
737 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
739 int cglo_flags = CGLO_GSTAT | CGLO_ENERGY;
740 compute_globals(gstat,
745 makeConstArrayRef(state->x),
746 makeConstArrayRef(state->v),
762 &observablesReducer);
763 // Clean up after pre-step use of compute_globals()
764 observablesReducer.markAsReadyToReduce();
767 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
768 the virial that should probably be addressed eventually. state->veta has better properies,
769 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
770 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
775 const bool bCalcEnerStep = true;
776 energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
782 ir->expandedvals.get(),
784 PTCouplingArrays({ state->boxv,
785 state->nosehoover_xi,
786 state->nosehoover_vxi,
788 state->nhpres_vxi }),
796 const bool do_ene = true;
797 const bool do_log = true;
799 const bool do_dr = ir->nstdisreout != 0;
800 const bool do_or = ir->nstorireout != 0;
802 EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
803 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
807 do_log ? fplog : nullptr,
815 pull_print_output(pull_work, step, t);
818 if (do_per_step(step, ir->nstlog))
820 if (fflush(fplog) != 0)
822 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
827 /* Print the remaining wall clock time for the run */
828 if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
832 fprintf(stderr, "\n");
834 print_time(stderr, walltime_accounting, step, ir, cr);
837 /* Ion/water position swapping.
838 * Not done in last step since trajectory writing happens before this call
839 * in the MD loop and exchanges would be lost anyway. */
840 if ((ir->eSwapCoords != SwapType::No) && (step > 0) && !isLastStep
841 && do_per_step(step, ir->swap->nstswap))
843 const bool doRerun = true;
852 MASTER(cr) && mdrunOptions.verbose,
858 /* read next frame from input trajectory */
859 isLastStep = !read_next_frame(oenv, status, &rerun_fr);
864 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
867 cycles = wallcycle_stop(wcycle, WallCycleCounter::Step);
868 if (haveDDAtomOrdering(*cr) && wcycle)
870 dd_cycles_add(cr->dd, cycles, ddCyclStep);
875 /* increase the MD step number */
879 observablesReducer.markAsReadyToReduce();
881 /* End of main MD loop */
883 /* Closing TNG files can include compressing data. Therefore it is good to do that
884 * before stopping the time measurements. */
885 mdoutf_tng_close(outf);
887 /* Stop measuring walltime */
888 walltime_accounting_end_time(walltime_accounting);
895 if (!thisRankHasDuty(cr, DUTY_PME))
897 /* Tell the PME only node to finish */
898 gmx_pme_send_finish(cr);
903 done_shellfc(fplog, shellfc, step_rel);
905 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);