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 "firstorderpressurecoupling.h"
82 #include "forceelement.h"
84 #include "nosehooverchains.h"
85 #include "parrinellorahmanbarostat.h"
86 #include "simulatoralgorithm.h"
87 #include "statepropagatordata.h"
88 #include "velocityscalingtemperaturecoupling.h"
92 void ModularSimulator::run()
94 GMX_LOG(legacySimulatorData_->mdlog.info)
96 .appendText("Using the modular simulator.");
98 ModularSimulatorAlgorithmBuilder algorithmBuilder(compat::make_not_null(legacySimulatorData_),
99 std::move(checkpointDataHolder_));
100 addIntegrationElements(&algorithmBuilder);
101 auto algorithm = algorithmBuilder.build();
103 while (const auto* task = algorithm.getNextTask())
110 void ModularSimulator::addIntegrationElements(ModularSimulatorAlgorithmBuilder* builder)
112 const bool isTrotter = inputrecNvtTrotter(legacySimulatorData_->inputrec)
113 || inputrecNptTrotter(legacySimulatorData_->inputrec)
114 || inputrecNphTrotter(legacySimulatorData_->inputrec);
115 if (legacySimulatorData_->inputrec->eI == IntegrationAlgorithm::MD)
117 // The leap frog integration algorithm
118 builder->add<ForceElement>();
119 builder->add<StatePropagatorData::Element>();
120 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::VRescale
121 || legacySimulatorData_->inputrec->etc == TemperatureCoupling::Berendsen
122 || legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
124 builder->add<VelocityScalingTemperatureCoupling>(Offset(-1),
126 ReportPreviousStepConservedEnergy::No,
127 PropagatorTag("LeapFrogPropagator"));
129 builder->add<Propagator<IntegrationStage::LeapFrog>>(
130 PropagatorTag("LeapFrogPropagator"), TimeStep(legacySimulatorData_->inputrec->delta_t));
131 if (legacySimulatorData_->constr)
133 builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
135 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>>();
136 if (legacySimulatorData_->inputrec->epc == PressureCoupling::ParrinelloRahman)
138 builder->add<ParrinelloRahmanBarostat>(Offset(-1), PropagatorTag("LeapFrogPropagator"));
140 else if (legacySimulatorData_->inputrec->epc == PressureCoupling::Berendsen
141 || legacySimulatorData_->inputrec->epc == PressureCoupling::CRescale)
143 builder->add<FirstOrderPressureCoupling>(0, ReportPreviousStepConservedEnergy::No);
146 else if (legacySimulatorData_->inputrec->eI == IntegrationAlgorithm::VV && !isTrotter)
148 // The velocity verlet integration algorithm
149 builder->add<ForceElement>();
150 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
151 PropagatorTag("VelocityHalfStep"), TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
152 if (legacySimulatorData_->constr)
154 builder->add<ConstraintsElement<ConstraintVariable::Velocities>>();
156 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
157 builder->add<StatePropagatorData::Element>();
158 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::VRescale
159 || legacySimulatorData_->inputrec->etc == TemperatureCoupling::Berendsen)
161 builder->add<VelocityScalingTemperatureCoupling>(
164 ReportPreviousStepConservedEnergy::Yes,
165 PropagatorTag("VelocityHalfAndPositionFullStep"));
167 else if (ETC_ANDERSEN(legacySimulatorData_->inputrec->etc))
169 builder->add<AndersenTemperatureCoupling>();
171 builder->add<Propagator<IntegrationStage::VelocityVerletPositionsAndVelocities>>(
172 PropagatorTag("VelocityHalfAndPositionFullStep"),
173 TimeStep(legacySimulatorData_->inputrec->delta_t));
174 if (legacySimulatorData_->constr)
176 builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
178 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
179 if (legacySimulatorData_->inputrec->epc == PressureCoupling::ParrinelloRahman)
181 builder->add<ParrinelloRahmanBarostat>(Offset(-1), PropagatorTag("VelocityHalfStep"));
183 else if (legacySimulatorData_->inputrec->epc == PressureCoupling::Berendsen
184 || legacySimulatorData_->inputrec->epc == PressureCoupling::CRescale)
186 builder->add<FirstOrderPressureCoupling>(0, ReportPreviousStepConservedEnergy::Yes);
189 else if (legacySimulatorData_->inputrec->eI == IntegrationAlgorithm::VV && isTrotter)
191 // For a new simulation, avoid the first Trotter half step
192 const auto scheduleTrotterFirstHalfOnInitStep =
193 ((legacySimulatorData_->startingBehavior == StartingBehavior::NewSimulation)
194 ? ScheduleOnInitStep::No
195 : ScheduleOnInitStep::Yes);
197 builder->add<ForceElement>();
198 // Propagate velocities from t-dt/2 to t
199 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
201 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
202 PropagatorTag("ScaleMTTKVPre1"));
204 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
205 PropagatorTag("VelocityHalfStep1"),
206 TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
207 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
209 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
210 PropagatorTag("ScaleMTTKVPost1"));
212 if (legacySimulatorData_->constr)
214 builder->add<ConstraintsElement<ConstraintVariable::Velocities>>();
216 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
218 // Propagate extended system variables from t-dt/2 to t
219 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
221 builder->add<MttkElement>(Offset(-1),
222 scheduleTrotterFirstHalfOnInitStep,
223 PropagatorTag("ScaleMTTKXPre"),
224 PropagatorTag("ScaleMTTKXPost"),
226 PropagatorTag("ScaleMTTKVPre1"),
227 PropagatorTag("ScaleMTTKVPost1"),
229 PropagatorTag("ScaleMTTKVPre2"),
230 PropagatorTag("ScaleMTTKVPost2"),
233 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
235 builder->add<NoseHooverChainsElement>(NhcUsage::System,
238 scheduleTrotterFirstHalfOnInitStep,
239 PropagatorTag("ScaleNHC"));
240 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(PropagatorTag("ScaleNHC"));
242 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
244 builder->add<NoseHooverChainsElement>(
245 NhcUsage::Barostat, Offset(-1), UseFullStepKE::Yes, scheduleTrotterFirstHalfOnInitStep);
247 // We have a full state at time t here
248 builder->add<StatePropagatorData::Element>();
250 // Propagate extended system variables from t to t+dt/2
251 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
253 builder->add<NoseHooverChainsElement>(
254 NhcUsage::Barostat, Offset(0), UseFullStepKE::Yes, ScheduleOnInitStep::Yes);
256 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
258 builder->add<NoseHooverChainsElement>(NhcUsage::System,
261 ScheduleOnInitStep::Yes,
262 PropagatorTag("VelocityHalfStep2"));
264 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
266 builder->add<MttkElement>(Offset(0),
267 ScheduleOnInitStep::Yes,
268 PropagatorTag("ScaleMTTKXPre"),
269 PropagatorTag("ScaleMTTKXPost"),
271 PropagatorTag("ScaleMTTKVPre1"),
272 PropagatorTag("ScaleMTTKVPost1"),
274 PropagatorTag("ScaleMTTKVPre2"),
275 PropagatorTag("ScaleMTTKVPost2"),
277 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
278 PropagatorTag("ScaleMTTKVPre2"));
281 // Propagate velocities from t to t+dt/2
282 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
283 PropagatorTag("VelocityHalfStep2"),
284 TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
285 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
287 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
288 PropagatorTag("ScaleMTTKVPost2"));
289 builder->add<Propagator<IntegrationStage::ScalePositions>>(
290 PropagatorTag("ScaleMTTKXPre"));
292 // Propagate positions from t to t+dt
293 builder->add<Propagator<IntegrationStage::PositionsOnly>>(
294 PropagatorTag("PositionFullStep"), TimeStep(legacySimulatorData_->inputrec->delta_t));
295 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
297 builder->add<Propagator<IntegrationStage::ScalePositions>>(
298 PropagatorTag("ScaleMTTKXPost"));
300 if (legacySimulatorData_->constr)
302 builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
304 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
306 // Propagate box from t to t+dt
307 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
309 builder->add<MttkBoxScaling>();
311 else if (legacySimulatorData_->inputrec->epc == PressureCoupling::CRescale)
313 // Legacy implementation allows combination of C-Rescale with Trotter Nose-Hoover
314 builder->add<FirstOrderPressureCoupling>(0, ReportPreviousStepConservedEnergy::Yes);
319 gmx_fatal(FARGS, "Integrator not implemented for the modular simulator.");
321 builder->add<EnergyData::Element>();
324 bool ModularSimulator::isInputCompatible(bool exitOnFailure,
325 const t_inputrec* inputrec,
327 const gmx_mtop_t& globalTopology,
328 const gmx_multisim_t* ms,
329 const ReplicaExchangeParameters& replExParams,
331 bool doEssentialDynamics,
334 auto conditionalAssert = [exitOnFailure](bool condition, const char* message) {
337 GMX_RELEASE_ASSERT(condition, message);
342 // GMX_USE_MODULAR_SIMULATOR allows to use modular simulator also for non-standard uses,
343 // such as the leap-frog integrator
344 const auto modularSimulatorExplicitlyTurnedOn = (getenv("GMX_USE_MODULAR_SIMULATOR") != nullptr);
345 // GMX_USE_MODULAR_SIMULATOR allows to use disable modular simulator for all uses,
346 // including the velocity-verlet integrator used by default
347 const auto modularSimulatorExplicitlyTurnedOff = (getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
350 !(modularSimulatorExplicitlyTurnedOn && modularSimulatorExplicitlyTurnedOff),
351 "Cannot have both GMX_USE_MODULAR_SIMULATOR=ON and GMX_DISABLE_MODULAR_SIMULATOR=ON. "
352 "Unset one of the two environment variables to explicitly chose which simulator to "
354 "or unset both to recover default behavior.");
357 !(modularSimulatorExplicitlyTurnedOff && inputrec->eI == IntegrationAlgorithm::VV
358 && inputrec->epc == PressureCoupling::ParrinelloRahman),
359 "Cannot use a Parrinello-Rahman barostat with md-vv and "
360 "GMX_DISABLE_MODULAR_SIMULATOR=ON, "
361 "as the Parrinello-Rahman barostat is not implemented in the legacy simulator. Unset "
362 "GMX_DISABLE_MODULAR_SIMULATOR or use a different pressure control algorithm.");
364 bool isInputCompatible = conditionalAssert(
365 inputrec->eI == IntegrationAlgorithm::MD || inputrec->eI == IntegrationAlgorithm::VV,
366 "Only integrators md and md-vv are supported by the modular simulator.");
367 isInputCompatible = isInputCompatible
368 && conditionalAssert(inputrec->eI != IntegrationAlgorithm::MD
369 || modularSimulatorExplicitlyTurnedOn,
370 "Set GMX_USE_MODULAR_SIMULATOR=ON to use the modular "
371 "simulator with integrator md.");
374 && conditionalAssert(
376 "Multiple time stepping is not supported by the modular simulator.");
379 && conditionalAssert(!doRerun, "Rerun is not supported by the modular simulator.");
382 && conditionalAssert(inputrec->efep == FreeEnergyPerturbationType::No
383 || inputrec->efep == FreeEnergyPerturbationType::Yes
384 || inputrec->efep == FreeEnergyPerturbationType::SlowGrowth,
385 "Expanded ensemble free energy calculation is not "
386 "supported by the modular simulator.");
387 isInputCompatible = isInputCompatible
388 && conditionalAssert(!inputrec->bPull,
389 "Pulling is not supported by the modular simulator.");
392 && conditionalAssert(inputrec->cos_accel == 0.0,
393 "Acceleration is not supported by the modular simulator.");
396 && conditionalAssert(!inputrecFrozenAtoms(inputrec),
397 "Freeze groups are not supported by the modular simulator.");
400 && conditionalAssert(
401 inputrec->deform[XX][XX] == 0.0 && inputrec->deform[XX][YY] == 0.0
402 && inputrec->deform[XX][ZZ] == 0.0 && inputrec->deform[YY][XX] == 0.0
403 && inputrec->deform[YY][YY] == 0.0 && inputrec->deform[YY][ZZ] == 0.0
404 && inputrec->deform[ZZ][XX] == 0.0 && inputrec->deform[ZZ][YY] == 0.0
405 && inputrec->deform[ZZ][ZZ] == 0.0,
406 "Deformation is not supported by the modular simulator.");
409 && conditionalAssert(gmx_mtop_interaction_count(globalTopology, IF_VSITE) == 0,
410 "Virtual sites are not supported by the modular simulator.");
411 isInputCompatible = isInputCompatible
412 && conditionalAssert(!inputrec->bDoAwh,
413 "AWH is not supported by the modular simulator.");
416 && conditionalAssert(gmx_mtop_ftype_count(globalTopology, F_DISRES) == 0,
417 "Distance restraints are not supported by the modular simulator.");
420 && conditionalAssert(
421 gmx_mtop_ftype_count(globalTopology, F_ORIRES) == 0,
422 "Orientation restraints are not supported by the modular simulator.");
425 && conditionalAssert(ms == nullptr,
426 "Multi-sim are not supported by the modular simulator.");
429 && conditionalAssert(replExParams.exchangeInterval == 0,
430 "Replica exchange is not supported by the modular simulator.");
432 int numEnsembleRestraintSystems;
435 numEnsembleRestraintSystems = fcd->disres->nsystems;
439 auto* distantRestraintEnsembleEnvVar = getenv("GMX_DISRE_ENSEMBLE_SIZE");
440 numEnsembleRestraintSystems =
441 (ms != nullptr && distantRestraintEnsembleEnvVar != nullptr)
442 ? static_cast<int>(strtol(distantRestraintEnsembleEnvVar, nullptr, 10))
447 && conditionalAssert(numEnsembleRestraintSystems <= 1,
448 "Ensemble restraints are not supported by the modular simulator.");
451 && conditionalAssert(!doSimulatedAnnealing(inputrec),
452 "Simulated annealing is not supported by the modular simulator.");
455 && conditionalAssert(!inputrec->bSimTemp,
456 "Simulated tempering is not supported by the modular simulator.");
457 isInputCompatible = isInputCompatible
458 && conditionalAssert(!inputrec->bExpanded,
459 "Expanded ensemble simulations are not supported by "
460 "the modular simulator.");
463 && conditionalAssert(!doEssentialDynamics,
464 "Essential dynamics is not supported by the modular simulator.");
465 isInputCompatible = isInputCompatible
466 && conditionalAssert(inputrec->eSwapCoords == SwapType::No,
467 "Ion / water position swapping is not supported by "
468 "the modular simulator.");
471 && conditionalAssert(!inputrec->bIMD,
472 "Interactive MD is not supported by the modular simulator.");
475 && conditionalAssert(!doMembed,
476 "Membrane embedding is not supported by the modular simulator.");
477 // TODO: Change this to the boolean passed when we merge the user interface change for the GPU update.
480 && conditionalAssert(
481 getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") == nullptr,
482 "Integration on the GPU is not supported by the modular simulator.");
483 // Modular simulator is centered around NS updates
484 // TODO: think how to handle nstlist == 0
485 isInputCompatible = isInputCompatible
486 && conditionalAssert(inputrec->nstlist != 0,
487 "Simulations without neighbor list update are not "
488 "supported by the modular simulator.");
489 isInputCompatible = isInputCompatible
490 && conditionalAssert(!GMX_FAHCORE,
491 "GMX_FAHCORE not supported by the modular simulator.");
492 if (!isInputCompatible
493 && (inputrec->eI == IntegrationAlgorithm::VV && inputrec->epc == PressureCoupling::ParrinelloRahman))
496 "Requested Parrinello-Rahman barostat with md-vv. This combination is only "
497 "available in the modular simulator. Some other selected options are, however, "
498 "only available in the legacy simulator. Use a different pressure control "
501 return isInputCompatible;
504 ModularSimulator::ModularSimulator(std::unique_ptr<LegacySimulatorData> legacySimulatorData,
505 std::unique_ptr<ReadCheckpointDataHolder> checkpointDataHolder) :
506 legacySimulatorData_(std::move(legacySimulatorData)),
507 checkpointDataHolder_(std::move(checkpointDataHolder))
509 checkInputForDisabledFunctionality();
512 void ModularSimulator::checkInputForDisabledFunctionality()
514 isInputCompatible(true,
515 legacySimulatorData_->inputrec,
516 legacySimulatorData_->mdrunOptions.rerun,
517 legacySimulatorData_->top_global,
518 legacySimulatorData_->ms,
519 legacySimulatorData_->replExParams,
520 legacySimulatorData_->fr->fcdata.get(),
521 opt2bSet("-ei", legacySimulatorData_->nfile, legacySimulatorData_->fnm),
522 legacySimulatorData_->membed != nullptr);
523 if (legacySimulatorData_->observablesHistory->edsamHistory)
526 "The checkpoint is from a run with essential dynamics sampling, "
527 "but the current run did not specify the -ei option. "
528 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
532 void ModularSimulator::readCheckpointToTrxFrame(t_trxframe* fr,
533 ReadCheckpointDataHolder* readCheckpointDataHolder,
534 const CheckpointHeaderContents& checkpointHeaderContents)
536 GMX_RELEASE_ASSERT(checkpointHeaderContents.isModularSimulatorCheckpoint,
537 "ModularSimulator::readCheckpointToTrxFrame can only read checkpoints "
538 "written by modular simulator.");
540 fr->step = int64_to_int(checkpointHeaderContents.step, "conversion of checkpoint to trajectory");
542 fr->time = checkpointHeaderContents.t;
546 StatePropagatorData::readCheckpointToTrxFrame(
547 fr, readCheckpointDataHolder->checkpointData(StatePropagatorData::checkpointID()));
548 if (readCheckpointDataHolder->keyExists(FreeEnergyPerturbationData::checkpointID()))
550 FreeEnergyPerturbationData::readCheckpointToTrxFrame(
551 fr, readCheckpointDataHolder->checkpointData(FreeEnergyPerturbationData::checkpointID()));
555 FreeEnergyPerturbationData::readCheckpointToTrxFrame(fr, std::nullopt);