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 "pullelement.h"
88 #include "simulatoralgorithm.h"
89 #include "statepropagatordata.h"
90 #include "velocityscalingtemperaturecoupling.h"
94 void ModularSimulator::run()
96 GMX_LOG(legacySimulatorData_->mdlog.info)
98 .appendText("Using the modular simulator.");
100 ModularSimulatorAlgorithmBuilder algorithmBuilder(compat::make_not_null(legacySimulatorData_),
101 std::move(checkpointDataHolder_));
102 addIntegrationElements(&algorithmBuilder);
103 auto algorithm = algorithmBuilder.build();
105 while (const auto* task = algorithm.getNextTask())
112 void ModularSimulator::addIntegrationElements(ModularSimulatorAlgorithmBuilder* builder)
114 const bool isTrotter = inputrecNvtTrotter(legacySimulatorData_->inputrec)
115 || inputrecNptTrotter(legacySimulatorData_->inputrec)
116 || inputrecNphTrotter(legacySimulatorData_->inputrec);
117 if (legacySimulatorData_->inputrec->eI == IntegrationAlgorithm::MD)
119 // The leap frog integration algorithm
120 builder->add<ForceElement>();
121 builder->add<StatePropagatorData::Element>();
122 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::VRescale
123 || legacySimulatorData_->inputrec->etc == TemperatureCoupling::Berendsen
124 || legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
126 builder->add<VelocityScalingTemperatureCoupling>(Offset(-1),
128 ReportPreviousStepConservedEnergy::No,
129 PropagatorTag("LeapFrogPropagator"));
131 builder->add<Propagator<IntegrationStage::LeapFrog>>(
132 PropagatorTag("LeapFrogPropagator"), TimeStep(legacySimulatorData_->inputrec->delta_t));
133 if (legacySimulatorData_->constr)
135 builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
138 if (legacySimulatorData_->inputrec->bPull)
140 builder->add<PullElement>();
143 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>>();
144 if (legacySimulatorData_->inputrec->epc == PressureCoupling::ParrinelloRahman)
146 builder->add<ParrinelloRahmanBarostat>(Offset(-1), PropagatorTag("LeapFrogPropagator"));
148 else if (legacySimulatorData_->inputrec->epc == PressureCoupling::Berendsen
149 || legacySimulatorData_->inputrec->epc == PressureCoupling::CRescale)
151 builder->add<FirstOrderPressureCoupling>(0, ReportPreviousStepConservedEnergy::No);
154 else if (legacySimulatorData_->inputrec->eI == IntegrationAlgorithm::VV && !isTrotter)
156 // The velocity verlet integration algorithm
157 builder->add<ForceElement>();
158 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
159 PropagatorTag("VelocityHalfStep"), TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
160 if (legacySimulatorData_->constr)
162 builder->add<ConstraintsElement<ConstraintVariable::Velocities>>();
164 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
165 // Here, we have x / v / f at the full time step
166 builder->add<StatePropagatorData::Element>();
167 if (legacySimulatorData_->inputrec->bExpanded)
169 builder->add<ExpandedEnsembleElement>();
171 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::VRescale
172 || legacySimulatorData_->inputrec->etc == TemperatureCoupling::Berendsen)
174 builder->add<VelocityScalingTemperatureCoupling>(
177 ReportPreviousStepConservedEnergy::Yes,
178 PropagatorTag("VelocityHalfAndPositionFullStep"));
180 else if (ETC_ANDERSEN(legacySimulatorData_->inputrec->etc))
182 builder->add<AndersenTemperatureCoupling>();
184 builder->add<Propagator<IntegrationStage::VelocityVerletPositionsAndVelocities>>(
185 PropagatorTag("VelocityHalfAndPositionFullStep"),
186 TimeStep(legacySimulatorData_->inputrec->delta_t));
187 if (legacySimulatorData_->constr)
189 builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
192 if (legacySimulatorData_->inputrec->bPull)
194 builder->add<PullElement>();
197 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
198 if (legacySimulatorData_->inputrec->epc == PressureCoupling::ParrinelloRahman)
200 builder->add<ParrinelloRahmanBarostat>(Offset(-1), PropagatorTag("VelocityHalfStep"));
202 else if (legacySimulatorData_->inputrec->epc == PressureCoupling::Berendsen
203 || legacySimulatorData_->inputrec->epc == PressureCoupling::CRescale)
205 builder->add<FirstOrderPressureCoupling>(0, ReportPreviousStepConservedEnergy::Yes);
208 else if (legacySimulatorData_->inputrec->eI == IntegrationAlgorithm::VV && isTrotter)
210 // For a new simulation, avoid the first Trotter half step
211 const auto scheduleTrotterFirstHalfOnInitStep =
212 ((legacySimulatorData_->startingBehavior == StartingBehavior::NewSimulation)
213 ? ScheduleOnInitStep::No
214 : ScheduleOnInitStep::Yes);
215 // Define the tags and offsets for MTTK pressure scaling
216 const MttkPropagatorConnectionDetails mttkPropagatorConnectionDetails = {
217 PropagatorTag("ScaleMTTKXPre"), PropagatorTag("ScaleMTTKXPost"), Offset(0),
218 PropagatorTag("ScaleMTTKVPre1"), PropagatorTag("ScaleMTTKVPost1"), Offset(1),
219 PropagatorTag("ScaleMTTKVPre2"), PropagatorTag("ScaleMTTKVPost2"), Offset(0)
222 builder->add<ForceElement>();
223 // Propagate velocities from t-dt/2 to t
224 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
226 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
227 PropagatorTag("ScaleMTTKVPre1"));
229 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
230 PropagatorTag("VelocityHalfStep1"),
231 TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
232 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
234 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
235 PropagatorTag("ScaleMTTKVPost1"));
237 if (legacySimulatorData_->constr)
239 builder->add<ConstraintsElement<ConstraintVariable::Velocities>>();
241 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
243 // Propagate extended system variables from t-dt/2 to t
244 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
246 builder->add<MttkElement>(
247 Offset(-1), scheduleTrotterFirstHalfOnInitStep, mttkPropagatorConnectionDetails);
249 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
251 builder->add<NoseHooverChainsElement>(NhcUsage::System,
254 scheduleTrotterFirstHalfOnInitStep,
255 PropagatorTag("ScaleNHC"));
256 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(PropagatorTag("ScaleNHC"));
258 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
260 builder->add<NoseHooverChainsElement>(NhcUsage::Barostat,
263 scheduleTrotterFirstHalfOnInitStep,
264 mttkPropagatorConnectionDetails);
266 // We have a full state at time t here
267 builder->add<StatePropagatorData::Element>();
268 if (legacySimulatorData_->inputrec->bExpanded)
270 builder->add<ExpandedEnsembleElement>();
273 // Propagate extended system variables from t to t+dt/2
274 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
276 builder->add<NoseHooverChainsElement>(NhcUsage::Barostat,
279 ScheduleOnInitStep::Yes,
280 mttkPropagatorConnectionDetails);
282 if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::NoseHoover)
284 builder->add<NoseHooverChainsElement>(NhcUsage::System,
287 ScheduleOnInitStep::Yes,
288 PropagatorTag("VelocityHalfStep2"));
290 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
292 builder->add<MttkElement>(Offset(0), ScheduleOnInitStep::Yes, mttkPropagatorConnectionDetails);
293 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
294 PropagatorTag("ScaleMTTKVPre2"));
297 // Propagate velocities from t to t+dt/2
298 builder->add<Propagator<IntegrationStage::VelocitiesOnly>>(
299 PropagatorTag("VelocityHalfStep2"),
300 TimeStep(0.5 * legacySimulatorData_->inputrec->delta_t));
301 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
303 builder->add<Propagator<IntegrationStage::ScaleVelocities>>(
304 PropagatorTag("ScaleMTTKVPost2"));
305 builder->add<Propagator<IntegrationStage::ScalePositions>>(
306 PropagatorTag("ScaleMTTKXPre"));
308 // Propagate positions from t to t+dt
309 builder->add<Propagator<IntegrationStage::PositionsOnly>>(
310 PropagatorTag("PositionFullStep"), TimeStep(legacySimulatorData_->inputrec->delta_t));
311 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
313 builder->add<Propagator<IntegrationStage::ScalePositions>>(
314 PropagatorTag("ScaleMTTKXPost"));
316 if (legacySimulatorData_->constr)
318 builder->add<ConstraintsElement<ConstraintVariable::Positions>>();
321 if (legacySimulatorData_->inputrec->bPull)
323 builder->add<PullElement>();
326 builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
328 // Propagate box from t to t+dt
329 if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
331 builder->add<MttkBoxScaling>(mttkPropagatorConnectionDetails);
333 else if (legacySimulatorData_->inputrec->epc == PressureCoupling::CRescale)
335 // Legacy implementation allows combination of C-Rescale with Trotter Nose-Hoover
336 builder->add<FirstOrderPressureCoupling>(0, ReportPreviousStepConservedEnergy::Yes);
341 gmx_fatal(FARGS, "Integrator not implemented for the modular simulator.");
343 builder->add<EnergyData::Element>();
346 bool ModularSimulator::isInputCompatible(bool exitOnFailure,
347 const t_inputrec* inputrec,
349 const gmx_mtop_t& globalTopology,
350 const gmx_multisim_t* ms,
351 const ReplicaExchangeParameters& replExParams,
353 bool doEssentialDynamics,
356 auto conditionalAssert = [exitOnFailure](bool condition, const char* message) {
359 GMX_RELEASE_ASSERT(condition, message);
364 // GMX_USE_MODULAR_SIMULATOR allows to use modular simulator also for non-standard uses,
365 // such as the leap-frog integrator
366 const auto modularSimulatorExplicitlyTurnedOn = (getenv("GMX_USE_MODULAR_SIMULATOR") != nullptr);
367 // GMX_USE_MODULAR_SIMULATOR allows to use disable modular simulator for all uses,
368 // including the velocity-verlet integrator used by default
369 const auto modularSimulatorExplicitlyTurnedOff = (getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
372 !(modularSimulatorExplicitlyTurnedOn && modularSimulatorExplicitlyTurnedOff),
373 "Cannot have both GMX_USE_MODULAR_SIMULATOR=ON and GMX_DISABLE_MODULAR_SIMULATOR=ON. "
374 "Unset one of the two environment variables to explicitly chose which simulator to "
376 "or unset both to recover default behavior.");
379 !(modularSimulatorExplicitlyTurnedOff && inputrec->eI == IntegrationAlgorithm::VV
380 && inputrec->epc == PressureCoupling::ParrinelloRahman),
381 "Cannot use a Parrinello-Rahman barostat with md-vv and "
382 "GMX_DISABLE_MODULAR_SIMULATOR=ON, "
383 "as the Parrinello-Rahman barostat is not implemented in the legacy simulator. Unset "
384 "GMX_DISABLE_MODULAR_SIMULATOR or use a different pressure control algorithm.");
386 bool isInputCompatible = conditionalAssert(
387 inputrec->eI == IntegrationAlgorithm::MD || inputrec->eI == IntegrationAlgorithm::VV,
388 "Only integrators md and md-vv are supported by the modular simulator.");
389 isInputCompatible = isInputCompatible
390 && conditionalAssert(inputrec->eI != IntegrationAlgorithm::MD
391 || modularSimulatorExplicitlyTurnedOn,
392 "Set GMX_USE_MODULAR_SIMULATOR=ON to use the modular "
393 "simulator with integrator md.");
396 && conditionalAssert(
398 "Multiple time stepping is not supported by the modular simulator.");
401 && conditionalAssert(!doRerun, "Rerun is not supported by the modular simulator.");
404 && conditionalAssert(!inputrec->useConstantAcceleration && inputrec->cos_accel == 0.0,
405 "Acceleration is not supported by the modular simulator.");
408 && conditionalAssert(!inputrecFrozenAtoms(inputrec),
409 "Freeze groups are not supported by the modular simulator.");
412 && conditionalAssert(
413 inputrec->deform[XX][XX] == 0.0 && inputrec->deform[XX][YY] == 0.0
414 && inputrec->deform[XX][ZZ] == 0.0 && inputrec->deform[YY][XX] == 0.0
415 && inputrec->deform[YY][YY] == 0.0 && inputrec->deform[YY][ZZ] == 0.0
416 && inputrec->deform[ZZ][XX] == 0.0 && inputrec->deform[ZZ][YY] == 0.0
417 && inputrec->deform[ZZ][ZZ] == 0.0,
418 "Deformation is not supported by the modular simulator.");
421 && conditionalAssert(gmx_mtop_interaction_count(globalTopology, IF_VSITE) == 0,
422 "Virtual sites are not supported by the modular simulator.");
423 isInputCompatible = isInputCompatible
424 && conditionalAssert(!inputrec->bDoAwh,
425 "AWH is not supported by the modular simulator.");
428 && conditionalAssert(gmx_mtop_ftype_count(globalTopology, F_DISRES) == 0,
429 "Distance restraints are not supported by the modular simulator.");
432 && conditionalAssert(
433 gmx_mtop_ftype_count(globalTopology, F_ORIRES) == 0,
434 "Orientation restraints are not supported by the modular simulator.");
437 && conditionalAssert(ms == nullptr,
438 "Multi-sim are not supported by the modular simulator.");
441 && conditionalAssert(replExParams.exchangeInterval == 0,
442 "Replica exchange is not supported by the modular simulator.");
444 int numEnsembleRestraintSystems;
447 numEnsembleRestraintSystems = fcd->disres->nsystems;
451 auto* distantRestraintEnsembleEnvVar = getenv("GMX_DISRE_ENSEMBLE_SIZE");
452 numEnsembleRestraintSystems =
453 (ms != nullptr && distantRestraintEnsembleEnvVar != nullptr)
454 ? static_cast<int>(strtol(distantRestraintEnsembleEnvVar, nullptr, 10))
459 && conditionalAssert(numEnsembleRestraintSystems <= 1,
460 "Ensemble restraints are not supported by the modular simulator.");
463 && conditionalAssert(!doSimulatedAnnealing(inputrec),
464 "Simulated annealing is not supported by the modular simulator.");
467 && conditionalAssert(!inputrec->bSimTemp,
468 "Simulated tempering is not supported by the modular simulator.");
471 && conditionalAssert(!doEssentialDynamics,
472 "Essential dynamics is not supported by the modular simulator.");
473 isInputCompatible = isInputCompatible
474 && conditionalAssert(inputrec->eSwapCoords == SwapType::No,
475 "Ion / water position swapping is not supported by "
476 "the modular simulator.");
479 && conditionalAssert(!inputrec->bIMD,
480 "Interactive MD is not supported by the modular simulator.");
483 && conditionalAssert(!doMembed,
484 "Membrane embedding is not supported by the modular simulator.");
485 // TODO: Change this to the boolean passed when we merge the user interface change for the GPU update.
488 && conditionalAssert(
489 getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") == nullptr,
490 "Integration on the GPU is not supported by the modular simulator.");
491 // Modular simulator is centered around NS updates
492 // TODO: think how to handle nstlist == 0
493 isInputCompatible = isInputCompatible
494 && conditionalAssert(inputrec->nstlist != 0,
495 "Simulations without neighbor list update are not "
496 "supported by the modular simulator.");
497 isInputCompatible = isInputCompatible
498 && conditionalAssert(!GMX_FAHCORE,
499 "GMX_FAHCORE not supported by the modular simulator.");
500 if (!isInputCompatible
501 && (inputrec->eI == IntegrationAlgorithm::VV && inputrec->epc == PressureCoupling::ParrinelloRahman))
504 "Requested Parrinello-Rahman barostat with md-vv. This combination is only "
505 "available in the modular simulator. Some other selected options are, however, "
506 "only available in the legacy simulator. Use a different pressure control "
509 return isInputCompatible;
512 ModularSimulator::ModularSimulator(std::unique_ptr<LegacySimulatorData> legacySimulatorData,
513 std::unique_ptr<ReadCheckpointDataHolder> checkpointDataHolder) :
514 legacySimulatorData_(std::move(legacySimulatorData)),
515 checkpointDataHolder_(std::move(checkpointDataHolder))
517 checkInputForDisabledFunctionality();
520 void ModularSimulator::checkInputForDisabledFunctionality()
522 isInputCompatible(true,
523 legacySimulatorData_->inputrec,
524 legacySimulatorData_->mdrunOptions.rerun,
525 legacySimulatorData_->top_global,
526 legacySimulatorData_->ms,
527 legacySimulatorData_->replExParams,
528 legacySimulatorData_->fr->fcdata.get(),
529 opt2bSet("-ei", legacySimulatorData_->nfile, legacySimulatorData_->fnm),
530 legacySimulatorData_->membed != nullptr);
531 if (legacySimulatorData_->observablesHistory->edsamHistory)
534 "The checkpoint is from a run with essential dynamics sampling, "
535 "but the current run did not specify the -ei option. "
536 "Either specify the -ei option to mdrun, or do not use this checkpoint file.");
540 void ModularSimulator::readCheckpointToTrxFrame(t_trxframe* fr,
541 ReadCheckpointDataHolder* readCheckpointDataHolder,
542 const CheckpointHeaderContents& checkpointHeaderContents)
544 GMX_RELEASE_ASSERT(checkpointHeaderContents.isModularSimulatorCheckpoint,
545 "ModularSimulator::readCheckpointToTrxFrame can only read checkpoints "
546 "written by modular simulator.");
548 fr->step = int64_to_int(checkpointHeaderContents.step, "conversion of checkpoint to trajectory");
550 fr->time = checkpointHeaderContents.t;
554 StatePropagatorData::readCheckpointToTrxFrame(
555 fr, readCheckpointDataHolder->checkpointData(StatePropagatorData::checkpointID()));
556 if (readCheckpointDataHolder->keyExists(FreeEnergyPerturbationData::checkpointID()))
558 FreeEnergyPerturbationData::readCheckpointToTrxFrame(
559 fr, readCheckpointDataHolder->checkpointData(FreeEnergyPerturbationData::checkpointID()));
563 FreeEnergyPerturbationData::readCheckpointToTrxFrame(fr, std::nullopt);