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/enerdata_utils.h"
50 #include "gromacs/mdlib/energyoutput.h"
51 #include "gromacs/mdlib/mdatoms.h"
52 #include "gromacs/mdlib/mdoutf.h"
53 #include "gromacs/mdlib/stat.h"
54 #include "gromacs/mdlib/tgroup.h"
55 #include "gromacs/mdlib/update.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/enerdata.h"
58 #include "gromacs/mdtypes/energyhistory.h"
59 #include "gromacs/mdtypes/group.h"
60 #include "gromacs/mdtypes/inputrec.h"
61 #include "gromacs/mdtypes/mdatom.h"
62 #include "gromacs/mdtypes/observableshistory.h"
63 #include "gromacs/mdtypes/pullhistory.h"
64 #include "gromacs/topology/topology.h"
66 #include "freeenergyperturbationdata.h"
67 #include "modularsimulator.h"
68 #include "simulatoralgorithm.h"
69 #include "statepropagatordata.h"
78 EnergyData::EnergyData(StatePropagatorData* statePropagatorData,
79 FreeEnergyPerturbationData* freeEnergyPerturbationData,
80 const gmx_mtop_t& globalTopology,
81 const t_inputrec* inputrec,
82 const MDAtoms* mdAtoms,
83 gmx_enerdata_t* enerd,
84 gmx_ekindata_t* ekind,
85 const Constraints* constr,
88 const MDModulesNotifiers& mdModulesNotifiers,
90 ObservablesHistory* observablesHistory,
91 StartingBehavior startingBehavior,
92 bool simulationsShareState) :
93 element_(std::make_unique<Element>(this, isMasterRank)),
94 isMasterRank_(isMasterRank),
99 needToSumEkinhOld_(false),
100 hasReadEkinFromCheckpoint_(false),
101 startingBehavior_(startingBehavior),
102 statePropagatorData_(statePropagatorData),
103 freeEnergyPerturbationData_(freeEnergyPerturbationData),
105 top_global_(globalTopology),
112 mdModulesNotifiers_(mdModulesNotifiers),
113 groups_(&globalTopology.groups),
114 observablesHistory_(observablesHistory),
115 simulationsShareState_(simulationsShareState)
117 clear_mat(forceVirial_);
118 clear_mat(shakeVirial_);
119 clear_mat(totalVirial_);
120 clear_mat(pressure_);
123 init_ekinstate(&ekinstate_, inputrec_);
124 observablesHistory_->energyHistory = std::make_unique<energyhistory_t>();
127 void EnergyData::Element::scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction)
133 auto writeEnergy = energyWritingStep_ == step;
134 auto isEnergyCalculationStep = energyCalculationStep_ == step;
135 auto isFreeEnergyCalculationStep = freeEnergyCalculationStep_ == step;
136 if (isEnergyCalculationStep || writeEnergy)
138 registerRunFunction([this, step, time, isEnergyCalculationStep, isFreeEnergyCalculationStep]() {
139 energyData_->doStep(step, time, isEnergyCalculationStep, isFreeEnergyCalculationStep);
144 registerRunFunction([this]() { energyData_->energyOutput_->recordNonEnergyStep(); });
148 void EnergyData::teardown()
150 if (inputrec_->nstcalcenergy > 0 && isMasterRank_)
152 energyOutput_->printEnergyConservation(fplog_, inputrec_->simulation_part, EI_MD(inputrec_->eI));
153 energyOutput_->printAverages(fplog_, groups_);
157 void EnergyData::Element::trajectoryWriterSetup(gmx_mdoutf* outf)
159 energyData_->setup(outf);
162 void EnergyData::setup(gmx_mdoutf* outf)
164 pull_t* pull_work = nullptr;
165 energyOutput_ = std::make_unique<EnergyOutput>(mdoutf_get_fp_ene(outf),
169 mdoutf_get_fp_dhdl(outf),
172 simulationsShareState_,
173 mdModulesNotifiers_);
180 initializeEnergyHistory(startingBehavior_, observablesHistory_, energyOutput_.get());
182 // TODO: This probably doesn't really belong here...
183 // but we have all we need in this element,
184 // so we'll leave it here for now!
185 double io = compute_io(inputrec_, top_global_.natoms, *groups_, energyOutput_->numEnergyTerms(), 1);
186 if ((io > 2000) && isMasterRank_)
188 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
190 if (!inputrec_->bContinuation)
192 real temp = enerd_->term[F_TEMP];
193 if (inputrec_->eI != IntegrationAlgorithm::VV)
195 /* Result of Ekin averaged over velocities of -half
196 * and +half step, while we only have -half step here.
200 fprintf(fplog_, "Initial temperature: %g K\n", temp);
204 std::optional<ITrajectoryWriterCallback> EnergyData::Element::registerTrajectoryWriterCallback(TrajectoryEvent event)
206 if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
208 return [this](gmx_mdoutf* mdoutf, Step step, Time time, bool writeTrajectory, bool writeLog) {
209 energyData_->write(mdoutf, step, time, writeTrajectory, writeLog);
215 std::optional<SignallerCallback> EnergyData::Element::registerTrajectorySignallerCallback(gmx::TrajectoryEvent event)
217 if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
219 return [this](Step step, Time /*unused*/) { energyWritingStep_ = step; };
224 std::optional<SignallerCallback> EnergyData::Element::registerEnergyCallback(EnergySignallerEvent event)
226 if (event == EnergySignallerEvent::EnergyCalculationStep && isMasterRank_)
228 return [this](Step step, Time /*unused*/) { energyCalculationStep_ = step; };
230 if (event == EnergySignallerEvent::FreeEnergyCalculationStep && isMasterRank_)
232 return [this](Step step, Time /*unused*/) { freeEnergyCalculationStep_ = step; };
237 void EnergyData::doStep(Step step, Time time, bool isEnergyCalculationStep, bool isFreeEnergyCalculationStep)
239 enerd_->term[F_ETOT] = enerd_->term[F_EPOT] + enerd_->term[F_EKIN];
240 if (freeEnergyPerturbationData_)
242 accumulateKineticLambdaComponents(
243 enerd_, freeEnergyPerturbationData_->constLambdaView(), *inputrec_->fepvals);
245 if (integratorHasConservedEnergyQuantity(inputrec_))
247 enerd_->term[F_ECONSERVED] = enerd_->term[F_ETOT];
248 for (const auto& energyContibution : conservedEnergyContributions_)
250 enerd_->term[F_ECONSERVED] += energyContibution(step, time);
253 matrix nullMatrix = {};
254 energyOutput_->addDataAtEnergyStep(
255 isFreeEnergyCalculationStep,
256 isEnergyCalculationStep,
258 mdAtoms_->mdatoms()->tmass,
260 inputrec_->fepvals.get(),
261 inputrec_->expandedvals.get(),
262 statePropagatorData_->constPreviousBox(),
263 PTCouplingArrays({ parrinelloRahmanBoxVelocities_ ? parrinelloRahmanBoxVelocities_() : nullMatrix,
268 freeEnergyPerturbationData_ ? freeEnergyPerturbationData_->currentFEPState() : 0,
278 void EnergyData::write(gmx_mdoutf* outf, Step step, Time time, bool writeTrajectory, bool writeLog)
282 energyOutput_->printHeader(fplog_, step, time);
285 bool do_dr = do_per_step(step, inputrec_->nstdisreout);
286 bool do_or = do_per_step(step, inputrec_->nstorireout);
288 // energyOutput_->printAnnealingTemperatures(writeLog ? fplog_ : nullptr, groups_, &(inputrec_->opts));
290 energyOutput_->printStepToEnergyFile(
291 mdoutf_get_fp_ene(outf), writeTrajectory, do_dr, do_or, writeLog ? fplog_ : nullptr, step, time, fcd_, awh);
294 void EnergyData::addToForceVirial(const tensor virial, Step step)
296 if (step > forceVirialStep_)
298 forceVirialStep_ = step;
299 clear_mat(forceVirial_);
301 m_add(forceVirial_, virial, forceVirial_);
304 void EnergyData::addToConstraintVirial(const tensor virial, Step step)
306 if (step > shakeVirialStep_)
308 shakeVirialStep_ = step;
309 clear_mat(shakeVirial_);
311 m_add(shakeVirial_, virial, shakeVirial_);
314 rvec* EnergyData::forceVirial(Step gmx_unused step)
316 if (step > forceVirialStep_)
318 forceVirialStep_ = step;
319 clear_mat(forceVirial_);
321 GMX_ASSERT(step >= forceVirialStep_ || forceVirialStep_ == -1,
322 "Asked for force virial of previous step.");
326 rvec* EnergyData::constraintVirial(Step gmx_unused step)
328 if (step > shakeVirialStep_)
330 shakeVirialStep_ = step;
331 clear_mat(shakeVirial_);
333 GMX_ASSERT(step >= shakeVirialStep_ || shakeVirialStep_ == -1,
334 "Asked for constraint virial of previous step.");
338 rvec* EnergyData::totalVirial(Step gmx_unused step)
340 if (step > totalVirialStep_)
342 totalVirialStep_ = step;
343 clear_mat(totalVirial_);
345 GMX_ASSERT(step >= totalVirialStep_ || totalVirialStep_ == -1,
346 "Asked for total virial of previous step.");
350 rvec* EnergyData::pressure(Step gmx_unused step)
352 if (step > pressureStep_)
354 pressureStep_ = step;
355 clear_mat(pressure_);
357 GMX_ASSERT(step >= pressureStep_ || pressureStep_ == -1, "Asked for pressure of previous step.");
361 real* EnergyData::muTot()
366 gmx_enerdata_t* EnergyData::enerdata()
371 gmx_ekindata_t* EnergyData::ekindata()
376 bool* EnergyData::needToSumEkinhOld()
378 return &needToSumEkinhOld_;
381 bool EnergyData::hasReadEkinFromCheckpoint() const
383 return hasReadEkinFromCheckpoint_;
389 * \brief Enum describing the contents EnergyData::Element writes to modular checkpoint
391 * When changing the checkpoint content, add a new element just above Count, and adjust the
392 * checkpoint functionality.
394 enum class CheckpointVersion
396 Base, //!< First version of modular checkpointing
397 Count //!< Number of entries. Add new versions right above this!
399 constexpr auto c_currentVersion = CheckpointVersion(int(CheckpointVersion::Count) - 1);
402 template<CheckpointDataOperation operation>
403 void EnergyData::Element::doCheckpointData(CheckpointData<operation>* checkpointData)
405 checkpointVersion(checkpointData, "EnergyData version", c_currentVersion);
407 energyData_->observablesHistory_->energyHistory->doCheckpoint<operation>(
408 checkpointData->subCheckpointData("energy history"));
409 energyData_->ekinstate_.doCheckpoint<operation>(checkpointData->subCheckpointData("ekinstate"));
412 void EnergyData::Element::saveCheckpointState(std::optional<WriteCheckpointData> checkpointData,
417 if (energyData_->needToSumEkinhOld_)
419 energyData_->ekinstate_.bUpToDate = false;
423 update_ekinstate(&energyData_->ekinstate_, energyData_->ekind_);
424 energyData_->ekinstate_.bUpToDate = true;
426 energyData_->energyOutput_->fillEnergyHistory(
427 energyData_->observablesHistory_->energyHistory.get());
428 doCheckpointData<CheckpointDataOperation::Write>(&checkpointData.value());
432 void EnergyData::Element::restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData,
437 doCheckpointData<CheckpointDataOperation::Read>(&checkpointData.value());
439 energyData_->hasReadEkinFromCheckpoint_ = MASTER(cr) ? energyData_->ekinstate_.bUpToDate : false;
442 gmx_bcast(sizeof(hasReadEkinFromCheckpoint_),
443 &energyData_->hasReadEkinFromCheckpoint_,
444 cr->mpi_comm_mygroup);
446 if (energyData_->hasReadEkinFromCheckpoint_)
448 // this takes care of broadcasting from master to agents
449 restore_ekinstate_from_state(cr, energyData_->ekind_, &energyData_->ekinstate_);
453 const std::string& EnergyData::Element::clientID()
458 void EnergyData::initializeEnergyHistory(StartingBehavior startingBehavior,
459 ObservablesHistory* observablesHistory,
460 EnergyOutput* energyOutput)
462 if (startingBehavior != StartingBehavior::NewSimulation)
464 /* Restore from energy history if appending to output files */
465 if (startingBehavior == StartingBehavior::RestartWithAppending)
467 /* If no history is available (because a checkpoint is from before
468 * it was written) make a new one later, otherwise restore it.
470 if (observablesHistory->energyHistory)
472 energyOutput->restoreFromEnergyHistory(*observablesHistory->energyHistory);
475 else if (observablesHistory->energyHistory)
477 /* We might have read an energy history from checkpoint.
478 * As we are not appending, we want to restart the statistics.
479 * Free the allocated memory and reset the counts.
481 observablesHistory->energyHistory = {};
482 /* We might have read a pull history from checkpoint.
483 * We will still want to keep the statistics, so that the files
484 * can be joined and still be meaningful.
485 * This means that observablesHistory_->pullHistory
486 * should not be reset.
490 if (!observablesHistory->energyHistory)
492 observablesHistory->energyHistory = std::make_unique<energyhistory_t>();
494 if (!observablesHistory->pullHistory)
496 observablesHistory->pullHistory = std::make_unique<PullHistory>();
498 /* Set the initial energy history */
499 energyOutput->fillEnergyHistory(observablesHistory->energyHistory.get());
502 void EnergyData::addConservedEnergyContribution(EnergyContribution&& energyContribution)
504 conservedEnergyContributions_.emplace_back(std::move(energyContribution));
507 void EnergyData::setParrinelloRahmanBoxVelocities(std::function<const rvec*()>&& parrinelloRahmanBoxVelocities)
509 GMX_RELEASE_ASSERT(!parrinelloRahmanBoxVelocities_,
510 "Received a second callback to the Parrinello-Rahman velocities");
511 parrinelloRahmanBoxVelocities_ = parrinelloRahmanBoxVelocities;
514 void EnergyData::updateKineticEnergy()
516 // The legacy sum_ekin function does not offer named types, so define variables for readability
517 // dEkin/dlambda is not handled here
518 real* dEkinDLambda = nullptr;
519 // Whether we use the full step kinetic energy (vs the average of half step KEs)
520 const bool useFullStepKineticEnergy = (inputrec_->eI == IntegrationAlgorithm::VV);
521 /* Whether we're ignoring the NHC scaling factor, only used if useFullStepKineticEnergy
522 * is true. (This parameter is confusing, as it is named `bScaleEkin`, but prompts the
523 * function to ignore scaling. There is no use case within modular simulator to ignore
524 * these, so we set this to false.) */
525 const bool ignoreScalingFactor = false;
527 enerd_->term[F_TEMP] = sum_ekin(
528 &(inputrec_->opts), ekind_, dEkinDLambda, useFullStepKineticEnergy, ignoreScalingFactor);
529 enerd_->term[F_EKIN] = trace(ekind_->ekin);
532 EnergyData::Element* EnergyData::element()
534 return element_.get();
537 EnergyData::Element::Element(EnergyData* energyData, bool isMasterRank) :
538 energyData_(energyData),
539 isMasterRank_(isMasterRank),
540 energyWritingStep_(-1),
541 energyCalculationStep_(-1),
542 freeEnergyCalculationStep_(-1)
546 ISimulatorElement* EnergyData::Element::getElementPointerImpl(
547 LegacySimulatorData gmx_unused* legacySimulatorData,
548 ModularSimulatorAlgorithmBuilderHelper gmx_unused* builderHelper,
549 StatePropagatorData gmx_unused* statePropagatorData,
550 EnergyData* energyData,
551 FreeEnergyPerturbationData gmx_unused* freeEnergyPerturbationData,
552 GlobalCommunicationHelper gmx_unused* globalCommunicationHelper)
554 return energyData->element();