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);
196 // Define the tags and offsets for MTTK pressure scaling
197 const MttkPropagatorConnectionDetails mttkPropagatorConnectionDetails = {
198 PropagatorTag("ScaleMTTKXPre"), PropagatorTag("ScaleMTTKXPost"), Offset(0),
199 PropagatorTag("ScaleMTTKVPre1"), PropagatorTag("ScaleMTTKVPost1"), Offset(1),
200 PropagatorTag("ScaleMTTKVPre2"), PropagatorTag("ScaleMTTKVPost2"), Offset(0)
203 builder->add<ForceElement>();
204 // Propagate velocities from t-dt/2 to t
205 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
207 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
208 PropagatorTag("ScaleMTTKVPre1"));
210 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
211 PropagatorTag("VelocityHalfStep1"),
212 TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
213 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
215 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
216 PropagatorTag("ScaleMTTKVPost1"));
218 if (legacySimulatorData_->constr)
220 builder->add<ConstraintsElement<ConstraintVariable::Velocities>>();
222 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
224 // Propagate extended system variables from t-dt/2 to t
225 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
227 builder->add<MttkElement>(
228 Offset(-1), scheduleTrotterFirstHalfOnInitStep, mttkPropagatorConnectionDetails);
230 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
232 builder->add<NoseHooverChainsElement>(NhcUsage::System,
235 scheduleTrotterFirstHalfOnInitStep,
236 PropagatorTag("ScaleNHC"));
237 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(PropagatorTag("ScaleNHC"));
239 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
241 builder->add<NoseHooverChainsElement>(NhcUsage::Barostat,
244 scheduleTrotterFirstHalfOnInitStep,
245 mttkPropagatorConnectionDetails);
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>(NhcUsage::Barostat,
256 ScheduleOnInitStep::Yes,
257 mttkPropagatorConnectionDetails);
259 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
261 builder->add<NoseHooverChainsElement>(NhcUsage::System,
264 ScheduleOnInitStep::Yes,
265 PropagatorTag("VelocityHalfStep2"));
267 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
269 builder->add<MttkElement>(Offset(0), ScheduleOnInitStep::Yes, mttkPropagatorConnectionDetails);
270 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
271 PropagatorTag("ScaleMTTKVPre2"));
274 // Propagate velocities from t to t+dt/2
275 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
276 PropagatorTag("VelocityHalfStep2"),
277 TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
278 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
280 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
281 PropagatorTag("ScaleMTTKVPost2"));
282 builder->add<Propagator<IntegrationStage::ScalePositions>>(
283 PropagatorTag("ScaleMTTKXPre"));
285 // Propagate positions from t to t+dt
286 builder->add<Propagator<IntegrationStage::PositionsOnly>>(
287 PropagatorTag("PositionFullStep"), TimeStep(legacySimulatorData_->inputrec->delta_t));
288 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
290 builder->add<Propagator<IntegrationStage::ScalePositions>>(
291 PropagatorTag("ScaleMTTKXPost"));
293 if (legacySimulatorData_->constr)
295 builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
297 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
299 // Propagate box from t to t+dt
300 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
302 builder->add<MttkBoxScaling>(mttkPropagatorConnectionDetails);
304 else if (legacySimulatorData_->inputrec->epc == PressureCoupling::CRescale)
306 // Legacy implementation allows combination of C-Rescale with Trotter Nose-Hoover
307 builder->add<FirstOrderPressureCoupling>(0, ReportPreviousStepConservedEnergy::Yes);
312 gmx_fatal(FARGS, "Integrator not implemented for the modular simulator.");
314 builder->add<EnergyData::Element>();
317 bool ModularSimulator::isInputCompatible(bool exitOnFailure,
318 const t_inputrec* inputrec,
320 const gmx_mtop_t& globalTopology,
321 const gmx_multisim_t* ms,
322 const ReplicaExchangeParameters& replExParams,
324 bool doEssentialDynamics,
327 auto conditionalAssert = [exitOnFailure](bool condition, const char* message) {
330 GMX_RELEASE_ASSERT(condition, message);
335 // GMX_USE_MODULAR_SIMULATOR allows to use modular simulator also for non-standard uses,
336 // such as the leap-frog integrator
337 const auto modularSimulatorExplicitlyTurnedOn = (getenv("GMX_USE_MODULAR_SIMULATOR") != nullptr);
338 // GMX_USE_MODULAR_SIMULATOR allows to use disable modular simulator for all uses,
339 // including the velocity-verlet integrator used by default
340 const auto modularSimulatorExplicitlyTurnedOff = (getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
343 !(modularSimulatorExplicitlyTurnedOn && modularSimulatorExplicitlyTurnedOff),
344 "Cannot have both GMX_USE_MODULAR_SIMULATOR=ON and GMX_DISABLE_MODULAR_SIMULATOR=ON. "
345 "Unset one of the two environment variables to explicitly chose which simulator to "
347 "or unset both to recover default behavior.");
350 !(modularSimulatorExplicitlyTurnedOff && inputrec->eI == IntegrationAlgorithm::VV
351 && inputrec->epc == PressureCoupling::ParrinelloRahman),
352 "Cannot use a Parrinello-Rahman barostat with md-vv and "
353 "GMX_DISABLE_MODULAR_SIMULATOR=ON, "
354 "as the Parrinello-Rahman barostat is not implemented in the legacy simulator. Unset "
355 "GMX_DISABLE_MODULAR_SIMULATOR or use a different pressure control algorithm.");
357 bool isInputCompatible = conditionalAssert(
358 inputrec->eI == IntegrationAlgorithm::MD || inputrec->eI == IntegrationAlgorithm::VV,
359 "Only integrators md and md-vv are supported by the modular simulator.");
360 isInputCompatible = isInputCompatible
361 && conditionalAssert(inputrec->eI != IntegrationAlgorithm::MD
362 || modularSimulatorExplicitlyTurnedOn,
363 "Set GMX_USE_MODULAR_SIMULATOR=ON to use the modular "
364 "simulator with integrator md.");
367 && conditionalAssert(
369 "Multiple time stepping is not supported by the modular simulator.");
372 && conditionalAssert(!doRerun, "Rerun is not supported by the modular simulator.");
375 && conditionalAssert(inputrec->efep == FreeEnergyPerturbationType::No
376 || inputrec->efep == FreeEnergyPerturbationType::Yes
377 || inputrec->efep == FreeEnergyPerturbationType::SlowGrowth,
378 "Expanded ensemble free energy calculation is not "
379 "supported by the modular simulator.");
380 isInputCompatible = isInputCompatible
381 && conditionalAssert(!inputrec->bPull,
382 "Pulling is not supported by the modular simulator.");
385 && conditionalAssert(inputrec->cos_accel == 0.0,
386 "Acceleration is not supported by the modular simulator.");
389 && conditionalAssert(!inputrecFrozenAtoms(inputrec),
390 "Freeze groups are not supported by the modular simulator.");
393 && conditionalAssert(
394 inputrec->deform[XX][XX] == 0.0 && inputrec->deform[XX][YY] == 0.0
395 && inputrec->deform[XX][ZZ] == 0.0 && inputrec->deform[YY][XX] == 0.0
396 && inputrec->deform[YY][YY] == 0.0 && inputrec->deform[YY][ZZ] == 0.0
397 && inputrec->deform[ZZ][XX] == 0.0 && inputrec->deform[ZZ][YY] == 0.0
398 && inputrec->deform[ZZ][ZZ] == 0.0,
399 "Deformation is not supported by the modular simulator.");
402 && conditionalAssert(gmx_mtop_interaction_count(globalTopology, IF_VSITE) == 0,
403 "Virtual sites are not supported by the modular simulator.");
404 isInputCompatible = isInputCompatible
405 && conditionalAssert(!inputrec->bDoAwh,
406 "AWH is not supported by the modular simulator.");
409 && conditionalAssert(gmx_mtop_ftype_count(globalTopology, F_DISRES) == 0,
410 "Distance restraints are not supported by the modular simulator.");
413 && conditionalAssert(
414 gmx_mtop_ftype_count(globalTopology, F_ORIRES) == 0,
415 "Orientation restraints are not supported by the modular simulator.");
418 && conditionalAssert(ms == nullptr,
419 "Multi-sim are not supported by the modular simulator.");
422 && conditionalAssert(replExParams.exchangeInterval == 0,
423 "Replica exchange is not supported by the modular simulator.");
425 int numEnsembleRestraintSystems;
428 numEnsembleRestraintSystems = fcd->disres->nsystems;
432 auto* distantRestraintEnsembleEnvVar = getenv("GMX_DISRE_ENSEMBLE_SIZE");
433 numEnsembleRestraintSystems =
434 (ms != nullptr && distantRestraintEnsembleEnvVar != nullptr)
435 ? static_cast<int>(strtol(distantRestraintEnsembleEnvVar, nullptr, 10))
440 && conditionalAssert(numEnsembleRestraintSystems <= 1,
441 "Ensemble restraints are not supported by the modular simulator.");
444 && conditionalAssert(!doSimulatedAnnealing(inputrec),
445 "Simulated annealing is not supported by the modular simulator.");
448 && conditionalAssert(!inputrec->bSimTemp,
449 "Simulated tempering is not supported by the modular simulator.");
450 isInputCompatible = isInputCompatible
451 && conditionalAssert(!inputrec->bExpanded,
452 "Expanded ensemble simulations are not supported by "
453 "the modular simulator.");
456 && conditionalAssert(!doEssentialDynamics,
457 "Essential dynamics is not supported by the modular simulator.");
458 isInputCompatible = isInputCompatible
459 && conditionalAssert(inputrec->eSwapCoords == SwapType::No,
460 "Ion / water position swapping is not supported by "
461 "the modular simulator.");
464 && conditionalAssert(!inputrec->bIMD,
465 "Interactive MD is not supported by the modular simulator.");
468 && conditionalAssert(!doMembed,
469 "Membrane embedding is not supported by the modular simulator.");
470 // TODO: Change this to the boolean passed when we merge the user interface change for the GPU update.
473 && conditionalAssert(
474 getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") == nullptr,
475 "Integration on the GPU is not supported by the modular simulator.");
476 // Modular simulator is centered around NS updates
477 // TODO: think how to handle nstlist == 0
478 isInputCompatible = isInputCompatible
479 && conditionalAssert(inputrec->nstlist != 0,
480 "Simulations without neighbor list update are not "
481 "supported by the modular simulator.");
482 isInputCompatible = isInputCompatible
483 && conditionalAssert(!GMX_FAHCORE,
484 "GMX_FAHCORE not supported by the modular simulator.");
485 if (!isInputCompatible
486 && (inputrec->eI == IntegrationAlgorithm::VV && inputrec->epc == PressureCoupling::ParrinelloRahman))
489 "Requested Parrinello-Rahman barostat with md-vv. This combination is only "
490 "available in the modular simulator. Some other selected options are, however, "
491 "only available in the legacy simulator. Use a different pressure control "
494 return isInputCompatible;
497 ModularSimulator::ModularSimulator(std::unique_ptr<LegacySimulatorData> legacySimulatorData,
498 std::unique_ptr<ReadCheckpointDataHolder> checkpointDataHolder) :
499 legacySimulatorData_(std::move(legacySimulatorData)),
500 checkpointDataHolder_(std::move(checkpointDataHolder))
502 checkInputForDisabledFunctionality();
505 void ModularSimulator::checkInputForDisabledFunctionality()
507 isInputCompatible(true,
508 legacySimulatorData_->inputrec,
509 legacySimulatorData_->mdrunOptions.rerun,
510 legacySimulatorData_->top_global,
511 legacySimulatorData_->ms,
512 legacySimulatorData_->replExParams,
513 legacySimulatorData_->fr->fcdata.get(),
514 opt2bSet("-ei", legacySimulatorData_->nfile, legacySimulatorData_->fnm),
515 legacySimulatorData_->membed != nullptr);
516 if (legacySimulatorData_->observablesHistory->edsamHistory)
519 "The checkpoint is from a run with essential dynamics sampling, "
520 "but the current run did not specify the -ei option. "
521 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
525 void ModularSimulator::readCheckpointToTrxFrame(t_trxframe* fr,
526 ReadCheckpointDataHolder* readCheckpointDataHolder,
527 const CheckpointHeaderContents& checkpointHeaderContents)
529 GMX_RELEASE_ASSERT(checkpointHeaderContents.isModularSimulatorCheckpoint,
530 "ModularSimulator::readCheckpointToTrxFrame can only read checkpoints "
531 "written by modular simulator.");
533 fr->step = int64_to_int(checkpointHeaderContents.step, "conversion of checkpoint to trajectory");
535 fr->time = checkpointHeaderContents.t;
539 StatePropagatorData::readCheckpointToTrxFrame(
540 fr, readCheckpointDataHolder->checkpointData(StatePropagatorData::checkpointID()));
541 if (readCheckpointDataHolder->keyExists(FreeEnergyPerturbationData::checkpointID()))
543 FreeEnergyPerturbationData::readCheckpointToTrxFrame(
544 fr, readCheckpointDataHolder->checkpointData(FreeEnergyPerturbationData::checkpointID()));
548 FreeEnergyPerturbationData::readCheckpointToTrxFrame(fr, std::nullopt);