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 /* Domain decomposition could incorrectly miss a bonded
159 interaction, but checking for that requires a global
160 communication stage, which does not otherwise happen in DD
161 code. So we do that alongside the first global energy reduction
162 after a new DD is made. These variables handle whether the
163 check happens, and the result it returns. */
164 bool shouldCheckNumberOfBondedInteractions = false;
165 int totalNumberOfBondedInteractions = -1;
167 SimulationSignals signals;
168 // Most global communnication stages don't propagate mdrun
169 // signals, and will use this object to achieve that.
170 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
174 gmx_fatal(FARGS, "Expanded ensemble not supported by MiMiC.");
178 gmx_fatal(FARGS, "Simulated tempering not supported by MiMiC.");
182 gmx_fatal(FARGS, "AWH not supported by MiMiC.");
184 if (replExParams.exchangeInterval > 0)
186 gmx_fatal(FARGS, "Replica exchange not supported by MiMiC.");
188 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
190 gmx_fatal(FARGS, "Essential dynamics not supported by MiMiC.");
194 gmx_fatal(FARGS, "Interactive MD not supported by MiMiC.");
198 gmx_fatal(FARGS, "Multiple simulations not supported by MiMiC.");
200 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc, [](SimulatedAnnealing i) {
201 return i != SimulatedAnnealing::No;
204 gmx_fatal(FARGS, "Simulated annealing not supported by MiMiC.");
207 /* Settings for rerun */
209 // TODO: Avoid changing inputrec (#3854)
210 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
211 nonConstInputrec->nstlist = 1;
212 nonConstInputrec->nstcalcenergy = 1;
213 nonConstInputrec->nstxout_compressed = 0;
215 int nstglobalcomm = 1;
216 const bool bNS = true;
220 MimicCommunicator::init();
221 auto nonConstGlobalTopology = const_cast<gmx_mtop_t*>(top_global);
222 MimicCommunicator::sendInitData(nonConstGlobalTopology, state_global->x);
223 // TODO: Avoid changing inputrec (#3854)
224 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
225 nonConstInputrec->nsteps = MimicCommunicator::getStepNumber();
227 if (DOMAINDECOMP(cr))
229 // TODO: Avoid changing inputrec (#3854)
230 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
231 gmx_bcast(sizeof(ir->nsteps), &nonConstInputrec->nsteps, cr->mpi_comm_mygroup);
234 const SimulationGroups* groups = &top_global->groups;
236 auto nonConstGlobalTopology = const_cast<gmx_mtop_t*>(top_global);
237 nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(*top_global);
240 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda);
242 const bool simulationsShareState = false;
243 gmx_mdoutf* outf = init_mdoutf(fplog,
254 StartingBehavior::NewSimulation,
255 simulationsShareState,
257 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
261 mdoutf_get_fp_dhdl(outf),
263 StartingBehavior::NewSimulation,
264 simulationsShareState,
267 gstat = global_stat_init(ir);
269 /* Check for polarizable models and flexible constraints */
270 shellfc = init_shell_flexcon(fplog,
272 constr ? constr->numFlexibleConstraints() : 0,
275 runScheduleWork->simulationWork.useGpuPme);
278 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
279 if ((io > 2000) && MASTER(cr))
281 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
285 // Local state only becomes valid now.
286 std::unique_ptr<t_state> stateInstance;
289 gmx_localtop_t top(top_global->ffparams);
291 if (DOMAINDECOMP(cr))
293 stateInstance = std::make_unique<t_state>();
294 state = stateInstance.get();
295 dd_init_local_state(*cr->dd, state_global, state);
297 /* Distribute the charge groups over the nodes from the master node */
298 dd_partition_system(fplog,
319 shouldCheckNumberOfBondedInteractions = true;
323 state_change_natoms(state_global, state_global->natoms);
324 /* Copy the pointer to the global state */
325 state = state_global;
327 mdAlgorithmsSetupAtomData(cr, *ir, *top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
330 auto mdatoms = mdAtoms->mdatoms();
332 // NOTE: The global state is no longer used at this point.
333 // But state_global is still used as temporary storage space for writing
334 // the global state to file and potentially for replica exchange.
335 // (Global topology should persist.)
337 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
339 if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
341 doFreeEnergyPerturbation = true;
347 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
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>{},
370 &totalNumberOfBondedInteractions,
374 checkNumberOfBondedInteractions(mdlog,
376 totalNumberOfBondedInteractions,
379 makeConstArrayRef(state->x),
381 &shouldCheckNumberOfBondedInteractions);
385 fprintf(stderr, "starting MiMiC MD run '%s'\n\n", *(top_global->name));
386 if (mdrunOptions.verbose)
389 "Calculated time to finish depends on nsteps from "
390 "run input file,\nwhich may not correspond to the time "
391 "needed to process input trajectory.\n\n");
393 fprintf(fplog, "\n");
396 walltime_accounting_start_time(walltime_accounting);
397 wallcycle_start(wcycle, ewcRUN);
398 print_start(fplog, cr, walltime_accounting, "mdrun");
400 /***********************************************************
404 ************************************************************/
411 "Simulations has constraints. Constraints will "
412 "be handled by CPMD.");
418 "MiMiC does not report kinetic energy, total energy, temperature, virial and "
421 step = ir->init_step;
424 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
425 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
429 mdrunOptions.reproducible,
431 mdrunOptions.maximumHoursToRun,
436 walltime_accounting);
438 // we don't do counter resetting in rerun - finish will always be valid
439 walltime_accounting_set_valid_finish(walltime_accounting);
441 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
443 /* and stop now if we should */
444 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
447 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel == ir->nsteps));
448 wallcycle_start(wcycle, ewcSTEP);
454 MimicCommunicator::getCoords(&state_global->x, state_global->natoms);
457 if (ir->efep != FreeEnergyPerturbationType::No)
459 state->lambda = currentLambdas(step, *(ir->fepvals), state_global->fep_state);
464 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
465 if (constructVsites && DOMAINDECOMP(cr))
468 "Vsite recalculation with -rerun is not implemented with domain "
470 "use a single rank");
474 wallcycle_start(wcycle, ewcVSITECONSTR);
475 vsite->construct(state->x, state->v, state->box, VSiteOperation::PositionsAndVelocities);
476 wallcycle_stop(wcycle, ewcVSITECONSTR);
480 if (DOMAINDECOMP(cr))
482 /* Repartition the domain decomposition */
483 const bool bMasterState = true;
484 dd_partition_system(fplog,
504 mdrunOptions.verbose);
505 shouldCheckNumberOfBondedInteractions = true;
510 EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
513 if (ir->efep != FreeEnergyPerturbationType::No)
515 update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
518 force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
519 | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 is solved
520 GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
524 /* Now is the time to relax the shells */
525 relax_shell_flexcon(fplog,
528 mdrunOptions.verbose,
540 state->x.arrayRefWithPadding(),
541 state->v.arrayRefWithPadding(),
556 ddBalanceRegionHandler);
560 /* The coordinates (x) are shifted (to get whole molecules)
562 * This is parallellized as well, and does communication too.
563 * Check comments in sim_util.c
566 gmx_edsam* ed = nullptr;
580 state->x.arrayRefWithPadding(),
593 GMX_FORCE_NS | force_flags,
594 ddBalanceRegionHandler);
597 /* Now we have the energies and forces corresponding to the
598 * coordinates at time t.
601 const bool isCheckpointingStep = false;
602 const bool doRerun = false;
603 const bool bSumEkinhOld = false;
604 do_md_trajectory_writing(fplog,
624 mdrunOptions.writeConfout,
628 stopHandler->setSignal();
631 const bool doInterSimSignal = false;
632 const bool doIntraSimSignal = true;
633 bool bSumEkinhOld = false;
634 t_vcm* vcm = nullptr;
635 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
637 compute_globals(gstat,
642 makeConstArrayRef(state->x),
643 makeConstArrayRef(state->v),
654 constr != nullptr ? constr->rmsdData() : gmx::ArrayRef<real>{},
657 &totalNumberOfBondedInteractions,
659 CGLO_GSTAT | CGLO_ENERGY
660 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
662 checkNumberOfBondedInteractions(mdlog,
664 totalNumberOfBondedInteractions,
667 makeConstArrayRef(state->x),
669 &shouldCheckNumberOfBondedInteractions);
673 gmx::HostVector<gmx::RVec> fglobal(top_global->natoms);
674 gmx::ArrayRef<gmx::RVec> ftemp;
675 gmx::ArrayRef<const gmx::RVec> flocal = f.view().force();
676 if (DOMAINDECOMP(cr))
678 ftemp = gmx::makeArrayRef(fglobal);
679 dd_collect_vec(cr->dd, state->ddp_count, state->ddp_count_cg_gl, state->cg_gl, flocal, ftemp);
683 ftemp = f.view().force();
688 MimicCommunicator::sendEnergies(enerd->term[F_EPOT]);
689 MimicCommunicator::sendForces(ftemp, state_global->natoms);
694 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
695 the virial that should probably be addressed eventually. state->veta has better properies,
696 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
697 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
699 if (ir->efep != FreeEnergyPerturbationType::No)
701 /* Sum up the foreign energy and dhdl terms for md and sd.
702 Currently done every step so that dhdl is correct in the .edr */
703 accumulateKineticLambdaComponents(enerd, state->lambda, *ir->fepvals);
709 const bool bCalcEnerStep = true;
710 energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
716 ir->expandedvals.get(),
718 PTCouplingArrays({ state->boxv,
719 state->nosehoover_xi,
720 state->nosehoover_vxi,
722 state->nhpres_vxi }),
732 const bool do_ene = true;
733 const bool do_log = true;
735 const bool do_dr = ir->nstdisreout != 0;
736 const bool do_or = ir->nstorireout != 0;
738 EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
739 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
743 do_log ? fplog : nullptr,
749 if (do_per_step(step, ir->nstlog))
751 if (fflush(fplog) != 0)
753 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
758 /* Print the remaining wall clock time for the run */
759 if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
763 fprintf(stderr, "\n");
765 print_time(stderr, walltime_accounting, step, ir, cr);
768 cycles = wallcycle_stop(wcycle, ewcSTEP);
769 if (DOMAINDECOMP(cr) && wcycle)
771 dd_cycles_add(cr->dd, cycles, ddCyclStep);
774 /* increase the MD step number */
778 /* End of main MD loop */
780 /* Closing TNG files can include compressing data. Therefore it is good to do that
781 * before stopping the time measurements. */
782 mdoutf_tng_close(outf);
784 /* Stop measuring walltime */
785 walltime_accounting_end_time(walltime_accounting);
789 MimicCommunicator::finalize();
792 if (!thisRankHasDuty(cr, DUTY_PME))
794 /* Tell the PME only node to finish */
795 gmx_pme_send_finish(cr);
800 done_shellfc(fplog, shellfc, step_rel);
802 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);