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.
38 * \brief Declares the loop for MiMiC QM/MM
40 * \author Viacheslav Bolnykh <v.bolnykh@hpc-leap.eu>
41 * \ingroup module_mdrun
53 #include "gromacs/applied_forces/awh/awh.h"
54 #include "gromacs/commandline/filenm.h"
55 #include "gromacs/domdec/collect.h"
56 #include "gromacs/domdec/dlbtiming.h"
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_network.h"
59 #include "gromacs/domdec/domdec_struct.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/enerdata.h"
107 #include "gromacs/mdtypes/energyhistory.h"
108 #include "gromacs/mdtypes/forcebuffers.h"
109 #include "gromacs/mdtypes/forcerec.h"
110 #include "gromacs/mdtypes/group.h"
111 #include "gromacs/mdtypes/inputrec.h"
112 #include "gromacs/mdtypes/interaction_const.h"
113 #include "gromacs/mdtypes/md_enums.h"
114 #include "gromacs/mdtypes/mdatom.h"
115 #include "gromacs/mdtypes/mdrunoptions.h"
116 #include "gromacs/mdtypes/observableshistory.h"
117 #include "gromacs/mdtypes/simulation_workload.h"
118 #include "gromacs/mdtypes/state.h"
119 #include "gromacs/mimic/communicator.h"
120 #include "gromacs/mimic/utilities.h"
121 #include "gromacs/pbcutil/pbc.h"
122 #include "gromacs/pulling/pull.h"
123 #include "gromacs/timing/wallcycle.h"
124 #include "gromacs/timing/walltime_accounting.h"
125 #include "gromacs/topology/atoms.h"
126 #include "gromacs/topology/idef.h"
127 #include "gromacs/topology/mtop_util.h"
128 #include "gromacs/topology/topology.h"
129 #include "gromacs/trajectory/trajectoryframe.h"
130 #include "gromacs/utility/basedefinitions.h"
131 #include "gromacs/utility/cstringutil.h"
132 #include "gromacs/utility/fatalerror.h"
133 #include "gromacs/utility/logger.h"
134 #include "gromacs/utility/real.h"
136 #include "legacysimulator.h"
137 #include "replicaexchange.h"
140 using gmx::SimulationSignaller;
142 void gmx::LegacySimulator::do_mimic()
144 const t_inputrec* ir = inputrec;
145 int64_t step, step_rel;
147 bool isLastStep = false;
148 bool doFreeEnergyPerturbation = false;
149 unsigned int force_flags;
150 tensor force_vir, shake_vir, total_vir, pres;
153 gmx_global_stat_t gstat;
154 gmx_shellfc_t* shellfc;
158 SimulationSignals signals;
159 // Most global communnication stages don't propagate mdrun
160 // signals, and will use this object to achieve that.
161 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
165 gmx_fatal(FARGS, "Expanded ensemble not supported by MiMiC.");
169 gmx_fatal(FARGS, "Simulated tempering not supported by MiMiC.");
173 gmx_fatal(FARGS, "AWH not supported by MiMiC.");
175 if (replExParams.exchangeInterval > 0)
177 gmx_fatal(FARGS, "Replica exchange not supported by MiMiC.");
179 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
181 gmx_fatal(FARGS, "Essential dynamics not supported by MiMiC.");
185 gmx_fatal(FARGS, "Interactive MD not supported by MiMiC.");
189 gmx_fatal(FARGS, "Multiple simulations not supported by MiMiC.");
191 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc, [](SimulatedAnnealing i) {
192 return i != SimulatedAnnealing::No;
195 gmx_fatal(FARGS, "Simulated annealing not supported by MiMiC.");
198 /* Settings for rerun */
200 // TODO: Avoid changing inputrec (#3854)
201 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
202 nonConstInputrec->nstlist = 1;
203 nonConstInputrec->nstcalcenergy = 1;
204 nonConstInputrec->nstxout_compressed = 0;
206 int nstglobalcomm = 1;
207 const bool bNS = true;
211 MimicCommunicator::init();
212 auto nonConstGlobalTopology = const_cast<gmx_mtop_t*>(&top_global);
213 MimicCommunicator::sendInitData(nonConstGlobalTopology, state_global->x);
214 // TODO: Avoid changing inputrec (#3854)
215 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
216 nonConstInputrec->nsteps = MimicCommunicator::getStepNumber();
218 if (DOMAINDECOMP(cr))
220 // TODO: Avoid changing inputrec (#3854)
221 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
222 gmx_bcast(sizeof(ir->nsteps), &nonConstInputrec->nsteps, cr->mpi_comm_mygroup);
225 const SimulationGroups* groups = &top_global.groups;
227 auto nonConstGlobalTopology = const_cast<gmx_mtop_t*>(&top_global);
228 nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(top_global);
231 initialize_lambdas(fplog,
235 ir->simtempvals->temperatures,
236 gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
238 &state_global->fep_state,
239 state_global->lambda);
241 const bool simulationsShareState = false;
242 gmx_mdoutf* outf = init_mdoutf(fplog,
253 StartingBehavior::NewSimulation,
254 simulationsShareState,
256 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
260 mdoutf_get_fp_dhdl(outf),
262 StartingBehavior::NewSimulation,
263 simulationsShareState,
266 gstat = global_stat_init(ir);
268 /* Check for polarizable models and flexible constraints */
269 shellfc = init_shell_flexcon(fplog,
271 constr ? constr->numFlexibleConstraints() : 0,
274 runScheduleWork->simulationWork.useGpuPme);
277 double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
278 if ((io > 2000) && MASTER(cr))
280 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
284 // Local state only becomes valid now.
285 std::unique_ptr<t_state> stateInstance;
288 gmx_localtop_t top(top_global.ffparams);
290 if (DOMAINDECOMP(cr))
292 stateInstance = std::make_unique<t_state>();
293 state = stateInstance.get();
294 dd_init_local_state(*cr->dd, state_global, state);
296 /* Distribute the charge groups over the nodes from the master node */
297 dd_partition_system(fplog,
321 state_change_natoms(state_global, state_global->natoms);
322 /* Copy the pointer to the global state */
323 state = state_global;
325 mdAlgorithmsSetupAtomData(cr, *ir, top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
328 auto mdatoms = mdAtoms->mdatoms();
330 // NOTE: The global state is no longer used at this point.
331 // But state_global is still used as temporary storage space for writing
332 // the global state to file and potentially for replica exchange.
333 // (Global topology should persist.)
335 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
337 if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
339 doFreeEnergyPerturbation = true;
343 int cglo_flags = CGLO_GSTAT;
344 if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
346 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
348 bool bSumEkinhOld = false;
349 t_vcm* vcm = nullptr;
350 compute_globals(gstat,
355 makeConstArrayRef(state->x),
356 makeConstArrayRef(state->v),
367 gmx::ArrayRef<real>{},
372 if (DOMAINDECOMP(cr))
374 checkNumberOfBondedInteractions(
375 mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
381 fprintf(stderr, "starting MiMiC MD run '%s'\n\n", *(top_global.name));
382 if (mdrunOptions.verbose)
385 "Calculated time to finish depends on nsteps from "
386 "run input file,\nwhich may not correspond to the time "
387 "needed to process input trajectory.\n\n");
389 fprintf(fplog, "\n");
392 walltime_accounting_start_time(walltime_accounting);
393 wallcycle_start(wcycle, WallCycleCounter::Run);
394 print_start(fplog, cr, walltime_accounting, "mdrun");
396 /***********************************************************
400 ************************************************************/
407 "Simulations has constraints. Constraints will "
408 "be handled by CPMD.");
414 "MiMiC does not report kinetic energy, total energy, temperature, virial and "
417 step = ir->init_step;
420 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
421 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
425 mdrunOptions.reproducible,
427 mdrunOptions.maximumHoursToRun,
432 walltime_accounting);
434 // we don't do counter resetting in rerun - finish will always be valid
435 walltime_accounting_set_valid_finish(walltime_accounting);
437 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
439 /* and stop now if we should */
440 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
443 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel == ir->nsteps));
444 wallcycle_start(wcycle, WallCycleCounter::Step);
450 MimicCommunicator::getCoords(&state_global->x, state_global->natoms);
453 if (ir->efep != FreeEnergyPerturbationType::No)
455 state->lambda = currentLambdas(step, *(ir->fepvals), state_global->fep_state);
460 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
461 if (constructVsites && DOMAINDECOMP(cr))
464 "Vsite recalculation with -rerun is not implemented with domain "
466 "use a single rank");
470 wallcycle_start(wcycle, WallCycleCounter::VsiteConstr);
471 vsite->construct(state->x, state->v, state->box, VSiteOperation::PositionsAndVelocities);
472 wallcycle_stop(wcycle, WallCycleCounter::VsiteConstr);
476 if (DOMAINDECOMP(cr))
478 /* Repartition the domain decomposition */
479 const bool bMasterState = true;
480 dd_partition_system(fplog,
500 mdrunOptions.verbose);
505 EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
508 if (ir->efep != FreeEnergyPerturbationType::No)
510 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
513 force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
514 | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 is solved
515 GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
519 /* Now is the time to relax the shells */
520 relax_shell_flexcon(fplog,
523 mdrunOptions.verbose,
535 state->x.arrayRefWithPadding(),
536 state->v.arrayRefWithPadding(),
551 ddBalanceRegionHandler);
555 /* The coordinates (x) are shifted (to get whole molecules)
557 * This is parallellized as well, and does communication too.
558 * Check comments in sim_util.c
561 gmx_edsam* ed = nullptr;
575 state->x.arrayRefWithPadding(),
588 GMX_FORCE_NS | force_flags,
589 ddBalanceRegionHandler);
592 /* Now we have the energies and forces corresponding to the
593 * coordinates at time t.
596 const bool isCheckpointingStep = false;
597 const bool doRerun = false;
598 const bool bSumEkinhOld = false;
599 do_md_trajectory_writing(fplog,
619 mdrunOptions.writeConfout,
623 stopHandler->setSignal();
626 const bool doInterSimSignal = false;
627 const bool doIntraSimSignal = true;
628 bool bSumEkinhOld = false;
629 t_vcm* vcm = nullptr;
630 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
632 int cglo_flags = CGLO_GSTAT | CGLO_ENERGY;
633 if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
635 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
637 compute_globals(gstat,
642 makeConstArrayRef(state->x),
643 makeConstArrayRef(state->v),
654 constr != nullptr ? constr->rmsdData() : gmx::ArrayRef<real>{},
659 if (DOMAINDECOMP(cr))
661 checkNumberOfBondedInteractions(
662 mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
667 gmx::HostVector<gmx::RVec> fglobal(top_global.natoms);
668 gmx::ArrayRef<gmx::RVec> ftemp;
669 gmx::ArrayRef<const gmx::RVec> flocal = f.view().force();
670 if (DOMAINDECOMP(cr))
672 ftemp = gmx::makeArrayRef(fglobal);
673 dd_collect_vec(cr->dd, state->ddp_count, state->ddp_count_cg_gl, state->cg_gl, flocal, ftemp);
677 ftemp = f.view().force();
682 MimicCommunicator::sendEnergies(enerd->term[F_EPOT]);
683 MimicCommunicator::sendForces(ftemp, state_global->natoms);
688 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
689 the virial that should probably be addressed eventually. state->veta has better properies,
690 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
691 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
693 if (ir->efep != FreeEnergyPerturbationType::No)
695 /* Sum up the foreign energy and dhdl terms for md and sd.
696 Currently done every step so that dhdl is correct in the .edr */
697 accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
703 const bool bCalcEnerStep = true;
704 energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
710 ir->expandedvals.get(),
712 PTCouplingArrays({ state->boxv,
713 state->nosehoover_xi,
714 state->nosehoover_vxi,
716 state->nhpres_vxi }),
726 const bool do_ene = true;
727 const bool do_log = true;
729 const bool do_dr = ir->nstdisreout != 0;
730 const bool do_or = ir->nstorireout != 0;
732 EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
733 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
737 do_log ? fplog : nullptr,
743 if (do_per_step(step, ir->nstlog))
745 if (fflush(fplog) != 0)
747 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
752 /* Print the remaining wall clock time for the run */
753 if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
757 fprintf(stderr, "\n");
759 print_time(stderr, walltime_accounting, step, ir, cr);
762 cycles = wallcycle_stop(wcycle, WallCycleCounter::Step);
763 if (DOMAINDECOMP(cr) && wcycle)
765 dd_cycles_add(cr->dd, cycles, ddCyclStep);
768 /* increase the MD step number */
772 /* End of main MD loop */
774 /* Closing TNG files can include compressing data. Therefore it is good to do that
775 * before stopping the time measurements. */
776 mdoutf_tng_close(outf);
778 /* Stop measuring walltime */
779 walltime_accounting_end_time(walltime_accounting);
783 MimicCommunicator::finalize();
786 if (!thisRankHasDuty(cr, DUTY_PME))
788 /* Tell the PME only node to finish */
789 gmx_pme_send_finish(cr);
794 done_shellfc(fplog, shellfc, step_rel);
796 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);