2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
36 * \brief Defines the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
44 #include "modularsimulator.h"
46 #include "gromacs/commandline/filenm.h"
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/ewald/pme.h"
49 #include "gromacs/ewald/pme_load_balancing.h"
50 #include "gromacs/ewald/pme_pp.h"
51 #include "gromacs/fileio/checkpoint.h"
52 #include "gromacs/gmxlib/nrnb.h"
53 #include "gromacs/listed_forces/listed_forces.h"
54 #include "gromacs/mdlib/checkpointhandler.h"
55 #include "gromacs/mdlib/constr.h"
56 #include "gromacs/mdlib/coupling.h"
57 #include "gromacs/mdlib/energyoutput.h"
58 #include "gromacs/mdlib/mdatoms.h"
59 #include "gromacs/mdlib/resethandler.h"
60 #include "gromacs/mdrun/replicaexchange.h"
61 #include "gromacs/mdrun/shellfc.h"
62 #include "gromacs/mdrunutility/handlerestart.h"
63 #include "gromacs/mdrunutility/printtime.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/fcdata.h"
66 #include "gromacs/mdtypes/forcerec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/mdatom.h"
69 #include "gromacs/mdtypes/mdrunoptions.h"
70 #include "gromacs/mdtypes/observableshistory.h"
71 #include "gromacs/nbnxm/nbnxm.h"
72 #include "gromacs/topology/mtop_util.h"
73 #include "gromacs/topology/topology.h"
74 #include "gromacs/trajectory/trajectoryframe.h"
75 #include "gromacs/utility/fatalerror.h"
76 #include "gromacs/utility/int64_to_int.h"
78 #include "andersentemperaturecoupling.h"
79 #include "computeglobalselement.h"
80 #include "constraintelement.h"
81 #include "expandedensembleelement.h"
82 #include "firstorderpressurecoupling.h"
83 #include "forceelement.h"
85 #include "nosehooverchains.h"
86 #include "parrinellorahmanbarostat.h"
87 #include "simulatoralgorithm.h"
88 #include "statepropagatordata.h"
89 #include "velocityscalingtemperaturecoupling.h"
93 void ModularSimulator::run()
95 GMX_LOG(legacySimulatorData_->mdlog.info)
97 .appendText("Using the modular simulator.");
99 ModularSimulatorAlgorithmBuilder algorithmBuilder(compat::make_not_null(legacySimulatorData_),
100 std::move(checkpointDataHolder_));
101 addIntegrationElements(&algorithmBuilder);
102 auto algorithm = algorithmBuilder.build();
104 while (const auto* task = algorithm.getNextTask())
111 void ModularSimulator::addIntegrationElements(ModularSimulatorAlgorithmBuilder* builder)
113 const bool isTrotter = inputrecNvtTrotter(legacySimulatorData_->inputrec)
114 || inputrecNptTrotter(legacySimulatorData_->inputrec)
115 || inputrecNphTrotter(legacySimulatorData_->inputrec);
116 if (legacySimulatorData_->inputrec->eI == IntegrationAlgorithm::MD)
118 // The leap frog integration algorithm
119 builder->add<ForceElement>();
120 builder->add<StatePropagatorData::Element>();
121 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::VRescale
122 || legacySimulatorData_->inputrec->etc == TemperatureCoupling::Berendsen
123 || legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
125 builder->add<VelocityScalingTemperatureCoupling>(Offset(-1),
127 ReportPreviousStepConservedEnergy::No,
128 PropagatorTag("LeapFrogPropagator"));
130 builder->add<Propagator<IntegrationStage::LeapFrog>>(
131 PropagatorTag("LeapFrogPropagator"), TimeStep(legacySimulatorData_->inputrec->delta_t));
132 if (legacySimulatorData_->constr)
134 builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
136 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>>();
137 if (legacySimulatorData_->inputrec->epc == PressureCoupling::ParrinelloRahman)
139 builder->add<ParrinelloRahmanBarostat>(Offset(-1), PropagatorTag("LeapFrogPropagator"));
141 else if (legacySimulatorData_->inputrec->epc == PressureCoupling::Berendsen
142 || legacySimulatorData_->inputrec->epc == PressureCoupling::CRescale)
144 builder->add<FirstOrderPressureCoupling>(0, ReportPreviousStepConservedEnergy::No);
147 else if (legacySimulatorData_->inputrec->eI == IntegrationAlgorithm::VV && !isTrotter)
149 // The velocity verlet integration algorithm
150 builder->add<ForceElement>();
151 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
152 PropagatorTag("VelocityHalfStep"), TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
153 if (legacySimulatorData_->constr)
155 builder->add<ConstraintsElement<ConstraintVariable::Velocities>>();
157 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
158 // Here, we have x / v / f at the full time step
159 builder->add<StatePropagatorData::Element>();
160 if (legacySimulatorData_->inputrec->bExpanded)
162 builder->add<ExpandedEnsembleElement>();
164 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::VRescale
165 || legacySimulatorData_->inputrec->etc == TemperatureCoupling::Berendsen)
167 builder->add<VelocityScalingTemperatureCoupling>(
170 ReportPreviousStepConservedEnergy::Yes,
171 PropagatorTag("VelocityHalfAndPositionFullStep"));
173 else if (ETC_ANDERSEN(legacySimulatorData_->inputrec->etc))
175 builder->add<AndersenTemperatureCoupling>();
177 builder->add<Propagator<IntegrationStage::VelocityVerletPositionsAndVelocities>>(
178 PropagatorTag("VelocityHalfAndPositionFullStep"),
179 TimeStep(legacySimulatorData_->inputrec->delta_t));
180 if (legacySimulatorData_->constr)
182 builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
184 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
185 if (legacySimulatorData_->inputrec->epc == PressureCoupling::ParrinelloRahman)
187 builder->add<ParrinelloRahmanBarostat>(Offset(-1), PropagatorTag("VelocityHalfStep"));
189 else if (legacySimulatorData_->inputrec->epc == PressureCoupling::Berendsen
190 || legacySimulatorData_->inputrec->epc == PressureCoupling::CRescale)
192 builder->add<FirstOrderPressureCoupling>(0, ReportPreviousStepConservedEnergy::Yes);
195 else if (legacySimulatorData_->inputrec->eI == IntegrationAlgorithm::VV && isTrotter)
197 // For a new simulation, avoid the first Trotter half step
198 const auto scheduleTrotterFirstHalfOnInitStep =
199 ((legacySimulatorData_->startingBehavior == StartingBehavior::NewSimulation)
200 ? ScheduleOnInitStep::No
201 : ScheduleOnInitStep::Yes);
202 // Define the tags and offsets for MTTK pressure scaling
203 const MttkPropagatorConnectionDetails mttkPropagatorConnectionDetails = {
204 PropagatorTag("ScaleMTTKXPre"), PropagatorTag("ScaleMTTKXPost"), Offset(0),
205 PropagatorTag("ScaleMTTKVPre1"), PropagatorTag("ScaleMTTKVPost1"), Offset(1),
206 PropagatorTag("ScaleMTTKVPre2"), PropagatorTag("ScaleMTTKVPost2"), Offset(0)
209 builder->add<ForceElement>();
210 // Propagate velocities from t-dt/2 to t
211 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
213 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
214 PropagatorTag("ScaleMTTKVPre1"));
216 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
217 PropagatorTag("VelocityHalfStep1"),
218 TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
219 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
221 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
222 PropagatorTag("ScaleMTTKVPost1"));
224 if (legacySimulatorData_->constr)
226 builder->add<ConstraintsElement<ConstraintVariable::Velocities>>();
228 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
230 // Propagate extended system variables from t-dt/2 to t
231 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
233 builder->add<MttkElement>(
234 Offset(-1), scheduleTrotterFirstHalfOnInitStep, mttkPropagatorConnectionDetails);
236 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
238 builder->add<NoseHooverChainsElement>(NhcUsage::System,
241 scheduleTrotterFirstHalfOnInitStep,
242 PropagatorTag("ScaleNHC"));
243 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(PropagatorTag("ScaleNHC"));
245 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
247 builder->add<NoseHooverChainsElement>(NhcUsage::Barostat,
250 scheduleTrotterFirstHalfOnInitStep,
251 mttkPropagatorConnectionDetails);
253 // We have a full state at time t here
254 builder->add<StatePropagatorData::Element>();
255 if (legacySimulatorData_->inputrec->bExpanded)
257 builder->add<ExpandedEnsembleElement>();
260 // Propagate extended system variables from t to t+dt/2
261 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
263 builder->add<NoseHooverChainsElement>(NhcUsage::Barostat,
266 ScheduleOnInitStep::Yes,
267 mttkPropagatorConnectionDetails);
269 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
271 builder->add<NoseHooverChainsElement>(NhcUsage::System,
274 ScheduleOnInitStep::Yes,
275 PropagatorTag("VelocityHalfStep2"));
277 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
279 builder->add<MttkElement>(Offset(0), ScheduleOnInitStep::Yes, mttkPropagatorConnectionDetails);
280 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
281 PropagatorTag("ScaleMTTKVPre2"));
284 // Propagate velocities from t to t+dt/2
285 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
286 PropagatorTag("VelocityHalfStep2"),
287 TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
288 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
290 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
291 PropagatorTag("ScaleMTTKVPost2"));
292 builder->add<Propagator<IntegrationStage::ScalePositions>>(
293 PropagatorTag("ScaleMTTKXPre"));
295 // Propagate positions from t to t+dt
296 builder->add<Propagator<IntegrationStage::PositionsOnly>>(
297 PropagatorTag("PositionFullStep"), TimeStep(legacySimulatorData_->inputrec->delta_t));
298 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
300 builder->add<Propagator<IntegrationStage::ScalePositions>>(
301 PropagatorTag("ScaleMTTKXPost"));
303 if (legacySimulatorData_->constr)
305 builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
307 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
309 // Propagate box from t to t+dt
310 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
312 builder->add<MttkBoxScaling>(mttkPropagatorConnectionDetails);
314 else if (legacySimulatorData_->inputrec->epc == PressureCoupling::CRescale)
316 // Legacy implementation allows combination of C-Rescale with Trotter Nose-Hoover
317 builder->add<FirstOrderPressureCoupling>(0, ReportPreviousStepConservedEnergy::Yes);
322 gmx_fatal(FARGS, "Integrator not implemented for the modular simulator.");
324 builder->add<EnergyData::Element>();
327 bool ModularSimulator::isInputCompatible(bool exitOnFailure,
328 const t_inputrec* inputrec,
330 const gmx_mtop_t& globalTopology,
331 const gmx_multisim_t* ms,
332 const ReplicaExchangeParameters& replExParams,
334 bool doEssentialDynamics,
337 auto conditionalAssert = [exitOnFailure](bool condition, const char* message) {
340 GMX_RELEASE_ASSERT(condition, message);
345 // GMX_USE_MODULAR_SIMULATOR allows to use modular simulator also for non-standard uses,
346 // such as the leap-frog integrator
347 const auto modularSimulatorExplicitlyTurnedOn = (getenv("GMX_USE_MODULAR_SIMULATOR") != nullptr);
348 // GMX_USE_MODULAR_SIMULATOR allows to use disable modular simulator for all uses,
349 // including the velocity-verlet integrator used by default
350 const auto modularSimulatorExplicitlyTurnedOff = (getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
353 !(modularSimulatorExplicitlyTurnedOn && modularSimulatorExplicitlyTurnedOff),
354 "Cannot have both GMX_USE_MODULAR_SIMULATOR=ON and GMX_DISABLE_MODULAR_SIMULATOR=ON. "
355 "Unset one of the two environment variables to explicitly chose which simulator to "
357 "or unset both to recover default behavior.");
360 !(modularSimulatorExplicitlyTurnedOff && inputrec->eI == IntegrationAlgorithm::VV
361 && inputrec->epc == PressureCoupling::ParrinelloRahman),
362 "Cannot use a Parrinello-Rahman barostat with md-vv and "
363 "GMX_DISABLE_MODULAR_SIMULATOR=ON, "
364 "as the Parrinello-Rahman barostat is not implemented in the legacy simulator. Unset "
365 "GMX_DISABLE_MODULAR_SIMULATOR or use a different pressure control algorithm.");
367 bool isInputCompatible = conditionalAssert(
368 inputrec->eI == IntegrationAlgorithm::MD || inputrec->eI == IntegrationAlgorithm::VV,
369 "Only integrators md and md-vv are supported by the modular simulator.");
370 isInputCompatible = isInputCompatible
371 && conditionalAssert(inputrec->eI != IntegrationAlgorithm::MD
372 || modularSimulatorExplicitlyTurnedOn,
373 "Set GMX_USE_MODULAR_SIMULATOR=ON to use the modular "
374 "simulator with integrator md.");
377 && conditionalAssert(
379 "Multiple time stepping is not supported by the modular simulator.");
382 && conditionalAssert(!doRerun, "Rerun is not supported by the modular simulator.");
383 isInputCompatible = isInputCompatible
384 && conditionalAssert(!inputrec->bPull,
385 "Pulling is not supported by the modular simulator.");
388 && conditionalAssert(inputrec->cos_accel == 0.0,
389 "Acceleration is not supported by the modular simulator.");
392 && conditionalAssert(!inputrecFrozenAtoms(inputrec),
393 "Freeze groups are not supported by the modular simulator.");
396 && conditionalAssert(
397 inputrec->deform[XX][XX] == 0.0 && inputrec->deform[XX][YY] == 0.0
398 && inputrec->deform[XX][ZZ] == 0.0 && inputrec->deform[YY][XX] == 0.0
399 && inputrec->deform[YY][YY] == 0.0 && inputrec->deform[YY][ZZ] == 0.0
400 && inputrec->deform[ZZ][XX] == 0.0 && inputrec->deform[ZZ][YY] == 0.0
401 && inputrec->deform[ZZ][ZZ] == 0.0,
402 "Deformation is not supported by the modular simulator.");
405 && conditionalAssert(gmx_mtop_interaction_count(globalTopology, IF_VSITE) == 0,
406 "Virtual sites are not supported by the modular simulator.");
407 isInputCompatible = isInputCompatible
408 && conditionalAssert(!inputrec->bDoAwh,
409 "AWH is not supported by the modular simulator.");
412 && conditionalAssert(gmx_mtop_ftype_count(globalTopology, F_DISRES) == 0,
413 "Distance restraints are not supported by the modular simulator.");
416 && conditionalAssert(
417 gmx_mtop_ftype_count(globalTopology, F_ORIRES) == 0,
418 "Orientation restraints are not supported by the modular simulator.");
421 && conditionalAssert(ms == nullptr,
422 "Multi-sim are not supported by the modular simulator.");
425 && conditionalAssert(replExParams.exchangeInterval == 0,
426 "Replica exchange is not supported by the modular simulator.");
428 int numEnsembleRestraintSystems;
431 numEnsembleRestraintSystems = fcd->disres->nsystems;
435 auto* distantRestraintEnsembleEnvVar = getenv("GMX_DISRE_ENSEMBLE_SIZE");
436 numEnsembleRestraintSystems =
437 (ms != nullptr && distantRestraintEnsembleEnvVar != nullptr)
438 ? static_cast<int>(strtol(distantRestraintEnsembleEnvVar, nullptr, 10))
443 && conditionalAssert(numEnsembleRestraintSystems <= 1,
444 "Ensemble restraints are not supported by the modular simulator.");
447 && conditionalAssert(!doSimulatedAnnealing(inputrec),
448 "Simulated annealing is not supported by the modular simulator.");
451 && conditionalAssert(!inputrec->bSimTemp,
452 "Simulated tempering is not supported by the modular simulator.");
455 && conditionalAssert(!doEssentialDynamics,
456 "Essential dynamics is not supported by the modular simulator.");
457 isInputCompatible = isInputCompatible
458 && conditionalAssert(inputrec->eSwapCoords == SwapType::No,
459 "Ion / water position swapping is not supported by "
460 "the modular simulator.");
463 && conditionalAssert(!inputrec->bIMD,
464 "Interactive MD is not supported by the modular simulator.");
467 && conditionalAssert(!doMembed,
468 "Membrane embedding is not supported by the modular simulator.");
469 // TODO: Change this to the boolean passed when we merge the user interface change for the GPU update.
472 && conditionalAssert(
473 getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") == nullptr,
474 "Integration on the GPU is not supported by the modular simulator.");
475 // Modular simulator is centered around NS updates
476 // TODO: think how to handle nstlist == 0
477 isInputCompatible = isInputCompatible
478 && conditionalAssert(inputrec->nstlist != 0,
479 "Simulations without neighbor list update are not "
480 "supported by the modular simulator.");
481 isInputCompatible = isInputCompatible
482 && conditionalAssert(!GMX_FAHCORE,
483 "GMX_FAHCORE not supported by the modular simulator.");
484 if (!isInputCompatible
485 && (inputrec->eI == IntegrationAlgorithm::VV && inputrec->epc == PressureCoupling::ParrinelloRahman))
488 "Requested Parrinello-Rahman barostat with md-vv. This combination is only "
489 "available in the modular simulator. Some other selected options are, however, "
490 "only available in the legacy simulator. Use a different pressure control "
493 return isInputCompatible;
496 ModularSimulator::ModularSimulator(std::unique_ptr<LegacySimulatorData> legacySimulatorData,
497 std::unique_ptr<ReadCheckpointDataHolder> checkpointDataHolder) :
498 legacySimulatorData_(std::move(legacySimulatorData)),
499 checkpointDataHolder_(std::move(checkpointDataHolder))
501 checkInputForDisabledFunctionality();
504 void ModularSimulator::checkInputForDisabledFunctionality()
506 isInputCompatible(true,
507 legacySimulatorData_->inputrec,
508 legacySimulatorData_->mdrunOptions.rerun,
509 legacySimulatorData_->top_global,
510 legacySimulatorData_->ms,
511 legacySimulatorData_->replExParams,
512 legacySimulatorData_->fr->fcdata.get(),
513 opt2bSet("-ei", legacySimulatorData_->nfile, legacySimulatorData_->fnm),
514 legacySimulatorData_->membed != nullptr);
515 if (legacySimulatorData_->observablesHistory->edsamHistory)
518 "The checkpoint is from a run with essential dynamics sampling, "
519 "but the current run did not specify the -ei option. "
520 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
524 void ModularSimulator::readCheckpointToTrxFrame(t_trxframe* fr,
525 ReadCheckpointDataHolder* readCheckpointDataHolder,
526 const CheckpointHeaderContents& checkpointHeaderContents)
528 GMX_RELEASE_ASSERT(checkpointHeaderContents.isModularSimulatorCheckpoint,
529 "ModularSimulator::readCheckpointToTrxFrame can only read checkpoints "
530 "written by modular simulator.");
532 fr->step = int64_to_int(checkpointHeaderContents.step, "conversion of checkpoint to trajectory");
534 fr->time = checkpointHeaderContents.t;
538 StatePropagatorData::readCheckpointToTrxFrame(
539 fr, readCheckpointDataHolder->checkpointData(StatePropagatorData::checkpointID()));
540 if (readCheckpointDataHolder->keyExists(FreeEnergyPerturbationData::checkpointID()))
542 FreeEnergyPerturbationData::readCheckpointToTrxFrame(
543 fr, readCheckpointDataHolder->checkpointData(FreeEnergyPerturbationData::checkpointID()));
547 FreeEnergyPerturbationData::readCheckpointToTrxFrame(fr, std::nullopt);