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 microstate for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
44 #include "energydata.h"
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdlib/compute_io.h"
49 #include "gromacs/mdlib/coupling.h"
50 #include "gromacs/mdlib/enerdata_utils.h"
51 #include "gromacs/mdlib/energyoutput.h"
52 #include "gromacs/mdlib/mdatoms.h"
53 #include "gromacs/mdlib/mdoutf.h"
54 #include "gromacs/mdlib/stat.h"
55 #include "gromacs/mdlib/update.h"
56 #include "gromacs/mdrunutility/handlerestart.h"
57 #include "gromacs/mdtypes/checkpointdata.h"
58 #include "gromacs/mdtypes/commrec.h"
59 #include "gromacs/mdtypes/enerdata.h"
60 #include "gromacs/mdtypes/energyhistory.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/mdatom.h"
63 #include "gromacs/mdtypes/observableshistory.h"
64 #include "gromacs/mdtypes/pullhistory.h"
65 #include "gromacs/topology/topology.h"
67 #include "freeenergyperturbationdata.h"
68 #include "modularsimulator.h"
69 #include "parrinellorahmanbarostat.h"
70 #include "simulatoralgorithm.h"
71 #include "statepropagatordata.h"
72 #include "velocityscalingtemperaturecoupling.h"
81 EnergyData::EnergyData(StatePropagatorData* statePropagatorData,
82 FreeEnergyPerturbationData* freeEnergyPerturbationData,
83 const gmx_mtop_t* globalTopology,
84 const t_inputrec* inputrec,
85 const MDAtoms* mdAtoms,
86 gmx_enerdata_t* enerd,
87 gmx_ekindata_t* ekind,
88 const Constraints* constr,
91 const MdModulesNotifier& mdModulesNotifier,
93 ObservablesHistory* observablesHistory,
94 StartingBehavior startingBehavior,
95 bool simulationsShareState) :
96 element_(std::make_unique<Element>(this, isMasterRank)),
97 isMasterRank_(isMasterRank),
100 totalVirialStep_(-1),
102 needToSumEkinhOld_(false),
103 hasReadEkinFromCheckpoint_(false),
104 startingBehavior_(startingBehavior),
105 statePropagatorData_(statePropagatorData),
106 freeEnergyPerturbationData_(freeEnergyPerturbationData),
107 velocityScalingTemperatureCoupling_(nullptr),
108 parrinelloRahmanBarostat_(nullptr),
110 top_global_(globalTopology),
117 mdModulesNotifier_(mdModulesNotifier),
118 groups_(&globalTopology->groups),
119 observablesHistory_(observablesHistory),
120 simulationsShareState_(simulationsShareState)
122 clear_mat(forceVirial_);
123 clear_mat(shakeVirial_);
124 clear_mat(totalVirial_);
125 clear_mat(pressure_);
128 init_ekinstate(&ekinstate_, inputrec_);
129 observablesHistory_->energyHistory = std::make_unique<energyhistory_t>();
132 void EnergyData::Element::scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction)
138 auto writeEnergy = energyWritingStep_ == step;
139 auto isEnergyCalculationStep = energyCalculationStep_ == step;
140 auto isFreeEnergyCalculationStep = freeEnergyCalculationStep_ == step;
141 if (isEnergyCalculationStep || writeEnergy)
143 registerRunFunction([this, time, isEnergyCalculationStep, isFreeEnergyCalculationStep]() {
144 energyData_->doStep(time, isEnergyCalculationStep, isFreeEnergyCalculationStep);
149 registerRunFunction([this]() { energyData_->energyOutput_->recordNonEnergyStep(); });
153 void EnergyData::teardown()
155 if (inputrec_->nstcalcenergy > 0 && isMasterRank_)
157 energyOutput_->printEnergyConservation(fplog_, inputrec_->simulation_part, EI_MD(inputrec_->eI));
158 energyOutput_->printAverages(fplog_, groups_);
162 void EnergyData::Element::trajectoryWriterSetup(gmx_mdoutf* outf)
164 energyData_->setup(outf);
167 void EnergyData::setup(gmx_mdoutf* outf)
169 pull_t* pull_work = nullptr;
170 energyOutput_ = std::make_unique<EnergyOutput>(mdoutf_get_fp_ene(outf),
174 mdoutf_get_fp_dhdl(outf),
177 simulationsShareState_,
185 initializeEnergyHistory(startingBehavior_, observablesHistory_, energyOutput_.get());
187 // TODO: This probably doesn't really belong here...
188 // but we have all we need in this element,
189 // so we'll leave it here for now!
190 double io = compute_io(inputrec_, top_global_->natoms, *groups_, energyOutput_->numEnergyTerms(), 1);
191 if ((io > 2000) && isMasterRank_)
193 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
195 if (!inputrec_->bContinuation)
197 real temp = enerd_->term[F_TEMP];
198 if (inputrec_->eI != IntegrationAlgorithm::VV)
200 /* Result of Ekin averaged over velocities of -half
201 * and +half step, while we only have -half step here.
205 fprintf(fplog_, "Initial temperature: %g K\n", temp);
209 std::optional<ITrajectoryWriterCallback> EnergyData::Element::registerTrajectoryWriterCallback(TrajectoryEvent event)
211 if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
213 return [this](gmx_mdoutf* mdoutf, Step step, Time time, bool writeTrajectory, bool writeLog) {
214 energyData_->write(mdoutf, step, time, writeTrajectory, writeLog);
220 std::optional<SignallerCallback> EnergyData::Element::registerTrajectorySignallerCallback(gmx::TrajectoryEvent event)
222 if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
224 return [this](Step step, Time /*unused*/) { energyWritingStep_ = step; };
229 std::optional<SignallerCallback> EnergyData::Element::registerEnergyCallback(EnergySignallerEvent event)
231 if (event == EnergySignallerEvent::EnergyCalculationStep && isMasterRank_)
233 return [this](Step step, Time /*unused*/) { energyCalculationStep_ = step; };
235 if (event == EnergySignallerEvent::FreeEnergyCalculationStep && isMasterRank_)
237 return [this](Step step, Time /*unused*/) { freeEnergyCalculationStep_ = step; };
242 void EnergyData::doStep(Time time, bool isEnergyCalculationStep, bool isFreeEnergyCalculationStep)
244 enerd_->term[F_ETOT] = enerd_->term[F_EPOT] + enerd_->term[F_EKIN];
245 if (freeEnergyPerturbationData_)
247 accumulateKineticLambdaComponents(
248 enerd_, freeEnergyPerturbationData_->constLambdaView(), *inputrec_->fepvals);
250 if (integratorHasConservedEnergyQuantity(inputrec_))
252 enerd_->term[F_ECONSERVED] =
254 + (velocityScalingTemperatureCoupling_
255 ? velocityScalingTemperatureCoupling_->conservedEnergyContribution()
257 + (parrinelloRahmanBarostat_ ? parrinelloRahmanBarostat_->conservedEnergyContribution() : 0);
259 matrix nullMatrix = {};
260 energyOutput_->addDataAtEnergyStep(
261 isFreeEnergyCalculationStep,
262 isEnergyCalculationStep,
264 mdAtoms_->mdatoms()->tmass,
266 inputrec_->fepvals.get(),
267 inputrec_->expandedvals.get(),
268 statePropagatorData_->constPreviousBox(),
269 PTCouplingArrays({ parrinelloRahmanBarostat_ ? parrinelloRahmanBarostat_->boxVelocities() : nullMatrix,
274 freeEnergyPerturbationData_ ? freeEnergyPerturbationData_->currentFEPState() : 0,
284 void EnergyData::write(gmx_mdoutf* outf, Step step, Time time, bool writeTrajectory, bool writeLog)
288 energyOutput_->printHeader(fplog_, step, time);
291 bool do_dr = do_per_step(step, inputrec_->nstdisreout);
292 bool do_or = do_per_step(step, inputrec_->nstorireout);
294 // energyOutput_->printAnnealingTemperatures(writeLog ? fplog_ : nullptr, groups_, &(inputrec_->opts));
296 energyOutput_->printStepToEnergyFile(
297 mdoutf_get_fp_ene(outf), writeTrajectory, do_dr, do_or, writeLog ? fplog_ : nullptr, step, time, fcd_, awh);
300 void EnergyData::addToForceVirial(const tensor virial, Step step)
302 if (step > forceVirialStep_)
304 forceVirialStep_ = step;
305 clear_mat(forceVirial_);
307 m_add(forceVirial_, virial, forceVirial_);
310 void EnergyData::addToConstraintVirial(const tensor virial, Step step)
312 if (step > shakeVirialStep_)
314 shakeVirialStep_ = step;
315 clear_mat(shakeVirial_);
317 m_add(shakeVirial_, virial, shakeVirial_);
320 rvec* EnergyData::forceVirial(Step gmx_unused step)
322 if (step > forceVirialStep_)
324 forceVirialStep_ = step;
325 clear_mat(forceVirial_);
327 GMX_ASSERT(step >= forceVirialStep_ || forceVirialStep_ == -1,
328 "Asked for force virial of previous step.");
332 rvec* EnergyData::constraintVirial(Step gmx_unused step)
334 if (step > shakeVirialStep_)
336 shakeVirialStep_ = step;
337 clear_mat(shakeVirial_);
339 GMX_ASSERT(step >= shakeVirialStep_ || shakeVirialStep_ == -1,
340 "Asked for constraint virial of previous step.");
344 rvec* EnergyData::totalVirial(Step gmx_unused step)
346 if (step > totalVirialStep_)
348 totalVirialStep_ = step;
349 clear_mat(totalVirial_);
351 GMX_ASSERT(step >= totalVirialStep_ || totalVirialStep_ == -1,
352 "Asked for total virial of previous step.");
356 rvec* EnergyData::pressure(Step gmx_unused step)
358 if (step > pressureStep_)
360 pressureStep_ = step;
361 clear_mat(pressure_);
363 GMX_ASSERT(step >= pressureStep_ || pressureStep_ == -1, "Asked for pressure of previous step.");
367 real* EnergyData::muTot()
372 gmx_enerdata_t* EnergyData::enerdata()
377 gmx_ekindata_t* EnergyData::ekindata()
382 bool* EnergyData::needToSumEkinhOld()
384 return &needToSumEkinhOld_;
387 bool EnergyData::hasReadEkinFromCheckpoint() const
389 return hasReadEkinFromCheckpoint_;
395 * \brief Enum describing the contents EnergyData::Element writes to modular checkpoint
397 * When changing the checkpoint content, add a new element just above Count, and adjust the
398 * checkpoint functionality.
400 enum class CheckpointVersion
402 Base, //!< First version of modular checkpointing
403 Count //!< Number of entries. Add new versions right above this!
405 constexpr auto c_currentVersion = CheckpointVersion(int(CheckpointVersion::Count) - 1);
408 template<CheckpointDataOperation operation>
409 void EnergyData::Element::doCheckpointData(CheckpointData<operation>* checkpointData)
411 checkpointVersion(checkpointData, "EnergyData version", c_currentVersion);
413 energyData_->observablesHistory_->energyHistory->doCheckpoint<operation>(
414 checkpointData->subCheckpointData("energy history"));
415 energyData_->ekinstate_.doCheckpoint<operation>(checkpointData->subCheckpointData("ekinstate"));
418 void EnergyData::Element::saveCheckpointState(std::optional<WriteCheckpointData> checkpointData,
423 if (energyData_->needToSumEkinhOld_)
425 energyData_->ekinstate_.bUpToDate = false;
429 update_ekinstate(&energyData_->ekinstate_, energyData_->ekind_);
430 energyData_->ekinstate_.bUpToDate = true;
432 energyData_->energyOutput_->fillEnergyHistory(
433 energyData_->observablesHistory_->energyHistory.get());
434 doCheckpointData<CheckpointDataOperation::Write>(&checkpointData.value());
438 void EnergyData::Element::restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData,
443 doCheckpointData<CheckpointDataOperation::Read>(&checkpointData.value());
445 energyData_->hasReadEkinFromCheckpoint_ = MASTER(cr) ? energyData_->ekinstate_.bUpToDate : false;
448 gmx_bcast(sizeof(hasReadEkinFromCheckpoint_),
449 &energyData_->hasReadEkinFromCheckpoint_,
450 cr->mpi_comm_mygroup);
452 if (energyData_->hasReadEkinFromCheckpoint_)
454 // this takes care of broadcasting from master to agents
455 restore_ekinstate_from_state(cr, energyData_->ekind_, &energyData_->ekinstate_);
459 const std::string& EnergyData::Element::clientID()
464 void EnergyData::initializeEnergyHistory(StartingBehavior startingBehavior,
465 ObservablesHistory* observablesHistory,
466 EnergyOutput* energyOutput)
468 if (startingBehavior != StartingBehavior::NewSimulation)
470 /* Restore from energy history if appending to output files */
471 if (startingBehavior == StartingBehavior::RestartWithAppending)
473 /* If no history is available (because a checkpoint is from before
474 * it was written) make a new one later, otherwise restore it.
476 if (observablesHistory->energyHistory)
478 energyOutput->restoreFromEnergyHistory(*observablesHistory->energyHistory);
481 else if (observablesHistory->energyHistory)
483 /* We might have read an energy history from checkpoint.
484 * As we are not appending, we want to restart the statistics.
485 * Free the allocated memory and reset the counts.
487 observablesHistory->energyHistory = {};
488 /* We might have read a pull history from checkpoint.
489 * We will still want to keep the statistics, so that the files
490 * can be joined and still be meaningful.
491 * This means that observablesHistory_->pullHistory
492 * should not be reset.
496 if (!observablesHistory->energyHistory)
498 observablesHistory->energyHistory = std::make_unique<energyhistory_t>();
500 if (!observablesHistory->pullHistory)
502 observablesHistory->pullHistory = std::make_unique<PullHistory>();
504 /* Set the initial energy history */
505 energyOutput->fillEnergyHistory(observablesHistory->energyHistory.get());
508 void EnergyData::setVelocityScalingTemperatureCoupling(const VelocityScalingTemperatureCoupling* velocityScalingTemperatureCoupling)
510 velocityScalingTemperatureCoupling_ = velocityScalingTemperatureCoupling;
513 void EnergyData::setParrinelloRahamnBarostat(const gmx::ParrinelloRahmanBarostat* parrinelloRahmanBarostat)
515 parrinelloRahmanBarostat_ = parrinelloRahmanBarostat;
518 EnergyData::Element* EnergyData::element()
520 return element_.get();
523 EnergyData::Element::Element(EnergyData* energyData, bool isMasterRank) :
524 energyData_(energyData),
525 isMasterRank_(isMasterRank),
526 energyWritingStep_(-1),
527 energyCalculationStep_(-1),
528 freeEnergyCalculationStep_(-1)
532 ISimulatorElement* EnergyData::Element::getElementPointerImpl(
533 LegacySimulatorData gmx_unused* legacySimulatorData,
534 ModularSimulatorAlgorithmBuilderHelper gmx_unused* builderHelper,
535 StatePropagatorData gmx_unused* statePropagatorData,
536 EnergyData* energyData,
537 FreeEnergyPerturbationData gmx_unused* freeEnergyPerturbationData,
538 GlobalCommunicationHelper gmx_unused* globalCommunicationHelper)
540 return energyData->element();