2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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/math/vec.h"
47 #include "gromacs/mdlib/compute_io.h"
48 #include "gromacs/mdlib/coupling.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/update.h"
55 #include "gromacs/mdrunutility/handlerestart.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/enerdata.h"
58 #include "gromacs/mdtypes/energyhistory.h"
59 #include "gromacs/mdtypes/inputrec.h"
60 #include "gromacs/mdtypes/mdatom.h"
61 #include "gromacs/mdtypes/observableshistory.h"
62 #include "gromacs/mdtypes/pullhistory.h"
63 #include "gromacs/mdtypes/state.h"
64 #include "gromacs/topology/topology.h"
66 #include "freeenergyperturbationdata.h"
67 #include "modularsimulator.h"
68 #include "parrinellorahmanbarostat.h"
69 #include "simulatoralgorithm.h"
70 #include "statepropagatordata.h"
71 #include "vrescalethermostat.h"
80 EnergyData::EnergyData(StatePropagatorData* statePropagatorData,
81 FreeEnergyPerturbationData* freeEnergyPerturbationData,
82 const gmx_mtop_t* globalTopology,
83 const t_inputrec* inputrec,
84 const MDAtoms* mdAtoms,
85 gmx_enerdata_t* enerd,
86 gmx_ekindata_t* ekind,
87 const Constraints* constr,
90 const MdModulesNotifier& mdModulesNotifier,
92 ObservablesHistory* observablesHistory,
93 StartingBehavior startingBehavior) :
94 element_(std::make_unique<Element>(this, isMasterRank)),
95 isMasterRank_(isMasterRank),
100 needToSumEkinhOld_(false),
101 startingBehavior_(startingBehavior),
102 statePropagatorData_(statePropagatorData),
103 freeEnergyPerturbationData_(freeEnergyPerturbationData),
104 vRescaleThermostat_(nullptr),
105 parrinelloRahmanBarostat_(nullptr),
107 top_global_(globalTopology),
114 mdModulesNotifier_(mdModulesNotifier),
115 groups_(&globalTopology->groups),
116 observablesHistory_(observablesHistory)
118 clear_mat(forceVirial_);
119 clear_mat(shakeVirial_);
120 clear_mat(totalVirial_);
121 clear_mat(pressure_);
124 if (freeEnergyPerturbationData_)
126 dummyLegacyState_.flags = (1U << estFEPSTATE);
130 void EnergyData::Element::scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction)
136 auto writeEnergy = energyWritingStep_ == step;
137 auto isEnergyCalculationStep = energyCalculationStep_ == step;
138 auto isFreeEnergyCalculationStep = freeEnergyCalculationStep_ == step;
139 if (isEnergyCalculationStep || writeEnergy)
141 registerRunFunction([this, time, isEnergyCalculationStep, isFreeEnergyCalculationStep]() {
142 energyData_->doStep(time, isEnergyCalculationStep, isFreeEnergyCalculationStep);
147 registerRunFunction([this]() { energyData_->energyOutput_->recordNonEnergyStep(); });
151 void EnergyData::teardown()
153 if (inputrec_->nstcalcenergy > 0 && isMasterRank_)
155 energyOutput_->printAverages(fplog_, groups_);
159 void EnergyData::Element::trajectoryWriterSetup(gmx_mdoutf* outf)
161 energyData_->setup(outf);
164 void EnergyData::setup(gmx_mdoutf* outf)
166 pull_t* pull_work = nullptr;
167 energyOutput_ = std::make_unique<EnergyOutput>(mdoutf_get_fp_ene(outf), top_global_, inputrec_,
168 pull_work, mdoutf_get_fp_dhdl(outf), false,
169 startingBehavior_, mdModulesNotifier_);
176 initializeEnergyHistory(startingBehavior_, observablesHistory_, energyOutput_.get());
178 // TODO: This probably doesn't really belong here...
179 // but we have all we need in this element,
180 // so we'll leave it here for now!
181 double io = compute_io(inputrec_, top_global_->natoms, *groups_, energyOutput_->numEnergyTerms(), 1);
182 if ((io > 2000) && isMasterRank_)
184 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
186 if (!inputrec_->bContinuation)
188 real temp = enerd_->term[F_TEMP];
189 if (inputrec_->eI != eiVV)
191 /* Result of Ekin averaged over velocities of -half
192 * and +half step, while we only have -half step here.
196 fprintf(fplog_, "Initial temperature: %g K\n", temp);
200 std::optional<ITrajectoryWriterCallback> EnergyData::Element::registerTrajectoryWriterCallback(TrajectoryEvent event)
202 if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
204 return [this](gmx_mdoutf* mdoutf, Step step, Time time, bool writeTrajectory, bool writeLog) {
205 energyData_->write(mdoutf, step, time, writeTrajectory, writeLog);
211 std::optional<SignallerCallback> EnergyData::Element::registerTrajectorySignallerCallback(gmx::TrajectoryEvent event)
213 if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
215 return [this](Step step, Time /*unused*/) { energyWritingStep_ = step; };
220 std::optional<SignallerCallback> EnergyData::Element::registerEnergyCallback(EnergySignallerEvent event)
222 if (event == EnergySignallerEvent::EnergyCalculationStep && isMasterRank_)
224 return [this](Step step, Time /*unused*/) { energyCalculationStep_ = step; };
226 if (event == EnergySignallerEvent::FreeEnergyCalculationStep && isMasterRank_)
228 return [this](Step step, Time /*unused*/) { freeEnergyCalculationStep_ = step; };
233 void EnergyData::doStep(Time time, bool isEnergyCalculationStep, bool isFreeEnergyCalculationStep)
235 enerd_->term[F_ETOT] = enerd_->term[F_EPOT] + enerd_->term[F_EKIN];
236 if (vRescaleThermostat_)
238 dummyLegacyState_.therm_integral = vRescaleThermostat_->thermostatIntegral();
240 if (freeEnergyPerturbationData_)
242 accumulateKineticLambdaComponents(enerd_, freeEnergyPerturbationData_->constLambdaView(),
243 *inputrec_->fepvals);
244 dummyLegacyState_.fep_state = freeEnergyPerturbationData_->currentFEPState();
246 if (parrinelloRahmanBarostat_)
248 copy_mat(parrinelloRahmanBarostat_->boxVelocities(), dummyLegacyState_.boxv);
249 copy_mat(statePropagatorData_->constBox(), dummyLegacyState_.box);
251 if (integratorHasConservedEnergyQuantity(inputrec_))
253 enerd_->term[F_ECONSERVED] =
254 enerd_->term[F_ETOT] + NPT_energy(inputrec_, &dummyLegacyState_, nullptr);
256 energyOutput_->addDataAtEnergyStep(isFreeEnergyCalculationStep, isEnergyCalculationStep, time,
257 mdAtoms_->mdatoms()->tmass, enerd_, &dummyLegacyState_,
258 inputrec_->fepvals, inputrec_->expandedvals,
259 statePropagatorData_->constPreviousBox(), shakeVirial_,
260 forceVirial_, totalVirial_, pressure_, ekind_, muTot_, constr_);
263 void EnergyData::write(gmx_mdoutf* outf, Step step, Time time, bool writeTrajectory, bool writeLog)
267 energyOutput_->printHeader(fplog_, step, time);
270 bool do_dr = do_per_step(step, inputrec_->nstdisreout);
271 bool do_or = do_per_step(step, inputrec_->nstorireout);
273 // energyOutput_->printAnnealingTemperatures(writeLog ? fplog_ : nullptr, groups_, &(inputrec_->opts));
275 energyOutput_->printStepToEnergyFile(mdoutf_get_fp_ene(outf), writeTrajectory, do_dr, do_or,
276 writeLog ? fplog_ : nullptr, step, time, fcd_, awh);
279 void EnergyData::addToForceVirial(const tensor virial, Step step)
281 if (step > forceVirialStep_)
283 forceVirialStep_ = step;
284 clear_mat(forceVirial_);
286 m_add(forceVirial_, virial, forceVirial_);
289 void EnergyData::addToConstraintVirial(const tensor virial, Step step)
291 if (step > shakeVirialStep_)
293 shakeVirialStep_ = step;
294 clear_mat(shakeVirial_);
296 m_add(shakeVirial_, virial, shakeVirial_);
299 rvec* EnergyData::forceVirial(Step gmx_unused step)
301 if (step > forceVirialStep_)
303 forceVirialStep_ = step;
304 clear_mat(forceVirial_);
306 GMX_ASSERT(step >= forceVirialStep_ || forceVirialStep_ == -1,
307 "Asked for force virial of previous step.");
311 rvec* EnergyData::constraintVirial(Step gmx_unused step)
313 if (step > shakeVirialStep_)
315 shakeVirialStep_ = step;
316 clear_mat(shakeVirial_);
318 GMX_ASSERT(step >= shakeVirialStep_ || shakeVirialStep_ == -1,
319 "Asked for constraint virial of previous step.");
323 rvec* EnergyData::totalVirial(Step gmx_unused step)
325 if (step > totalVirialStep_)
327 totalVirialStep_ = step;
328 clear_mat(totalVirial_);
330 GMX_ASSERT(step >= totalVirialStep_ || totalVirialStep_ == -1,
331 "Asked for total virial of previous step.");
335 rvec* EnergyData::pressure(Step gmx_unused step)
337 if (step > pressureStep_)
339 pressureStep_ = step;
340 clear_mat(pressure_);
342 GMX_ASSERT(step >= pressureStep_ || pressureStep_ == -1,
343 "Asked for pressure of previous step.");
347 real* EnergyData::muTot()
352 gmx_enerdata_t* EnergyData::enerdata()
357 gmx_ekindata_t* EnergyData::ekindata()
362 bool* EnergyData::needToSumEkinhOld()
364 return &needToSumEkinhOld_;
367 void EnergyData::Element::writeCheckpoint(t_state gmx_unused* localState, t_state* globalState)
371 if (energyData_->needToSumEkinhOld_)
373 globalState->ekinstate.bUpToDate = false;
377 update_ekinstate(&globalState->ekinstate, energyData_->ekind_);
378 globalState->ekinstate.bUpToDate = true;
380 energyData_->energyOutput_->fillEnergyHistory(
381 energyData_->observablesHistory_->energyHistory.get());
385 void EnergyData::initializeEnergyHistory(StartingBehavior startingBehavior,
386 ObservablesHistory* observablesHistory,
387 EnergyOutput* energyOutput)
389 if (startingBehavior != StartingBehavior::NewSimulation)
391 /* Restore from energy history if appending to output files */
392 if (startingBehavior == StartingBehavior::RestartWithAppending)
394 /* If no history is available (because a checkpoint is from before
395 * it was written) make a new one later, otherwise restore it.
397 if (observablesHistory->energyHistory)
399 energyOutput->restoreFromEnergyHistory(*observablesHistory->energyHistory);
402 else if (observablesHistory->energyHistory)
404 /* We might have read an energy history from checkpoint.
405 * As we are not appending, we want to restart the statistics.
406 * Free the allocated memory and reset the counts.
408 observablesHistory->energyHistory = {};
409 /* We might have read a pull history from checkpoint.
410 * We will still want to keep the statistics, so that the files
411 * can be joined and still be meaningful.
412 * This means that observablesHistory_->pullHistory
413 * should not be reset.
417 if (!observablesHistory->energyHistory)
419 observablesHistory->energyHistory = std::make_unique<energyhistory_t>();
421 if (!observablesHistory->pullHistory)
423 observablesHistory->pullHistory = std::make_unique<PullHistory>();
425 /* Set the initial energy history */
426 energyOutput->fillEnergyHistory(observablesHistory->energyHistory.get());
429 void EnergyData::setVRescaleThermostat(const gmx::VRescaleThermostat* vRescaleThermostat)
431 vRescaleThermostat_ = vRescaleThermostat;
432 if (vRescaleThermostat_)
434 dummyLegacyState_.flags |= (1U << estTHERM_INT);
438 void EnergyData::setParrinelloRahamnBarostat(const gmx::ParrinelloRahmanBarostat* parrinelloRahmanBarostat)
440 parrinelloRahmanBarostat_ = parrinelloRahmanBarostat;
441 if (parrinelloRahmanBarostat_)
443 dummyLegacyState_.flags |= (1U << estBOX) | (1U << estBOXV);
447 EnergyData::Element* EnergyData::element()
449 return element_.get();
452 EnergyData::Element::Element(EnergyData* energyData, bool isMasterRank) :
453 energyData_(energyData),
454 isMasterRank_(isMasterRank),
455 energyWritingStep_(-1),
456 energyCalculationStep_(-1),
457 freeEnergyCalculationStep_(-1)
461 ISimulatorElement* EnergyData::Element::getElementPointerImpl(
462 LegacySimulatorData gmx_unused* legacySimulatorData,
463 ModularSimulatorAlgorithmBuilderHelper gmx_unused* builderHelper,
464 StatePropagatorData gmx_unused* statePropagatorData,
465 EnergyData* energyData,
466 FreeEnergyPerturbationData gmx_unused* freeEnergyPerturbationData,
467 GlobalCommunicationHelper gmx_unused* globalCommunicationHelper)
469 return energyData->element();