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 "energyelement.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/enerdata.h"
57 #include "gromacs/mdtypes/energyhistory.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/mdtypes/mdatom.h"
60 #include "gromacs/mdtypes/observableshistory.h"
61 #include "gromacs/mdtypes/pullhistory.h"
62 #include "gromacs/mdtypes/state.h"
63 #include "gromacs/topology/topology.h"
65 #include "freeenergyperturbationelement.h"
66 #include "parrinellorahmanbarostat.h"
67 #include "statepropagatordata.h"
68 #include "vrescalethermostat.h"
77 EnergyElement::EnergyElement(StatePropagatorData* statePropagatorData,
78 FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
79 const gmx_mtop_t* globalTopology,
80 const t_inputrec* inputrec,
81 const MDAtoms* mdAtoms,
82 gmx_enerdata_t* enerd,
83 gmx_ekindata_t* ekind,
84 const Constraints* constr,
87 const MdModulesNotifier& mdModulesNotifier,
89 ObservablesHistory* observablesHistory,
90 StartingBehavior startingBehavior) :
91 isMasterRank_(isMasterRank),
92 energyWritingStep_(-1),
93 energyCalculationStep_(-1),
94 freeEnergyCalculationStep_(-1),
99 needToSumEkinhOld_(false),
100 startingBehavior_(startingBehavior),
101 statePropagatorData_(statePropagatorData),
102 freeEnergyPerturbationElement_(freeEnergyPerturbationElement),
103 vRescaleThermostat_(nullptr),
104 parrinelloRahmanBarostat_(nullptr),
106 top_global_(globalTopology),
113 mdModulesNotifier_(mdModulesNotifier),
114 groups_(&globalTopology->groups),
115 observablesHistory_(observablesHistory)
117 clear_mat(forceVirial_);
118 clear_mat(shakeVirial_);
119 clear_mat(totalVirial_);
120 clear_mat(pressure_);
123 if (freeEnergyPerturbationElement_)
125 dummyLegacyState_.flags = (1U << estFEPSTATE);
129 void EnergyElement::scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction)
135 auto writeEnergy = energyWritingStep_ == step;
136 auto isEnergyCalculationStep = energyCalculationStep_ == step;
137 auto isFreeEnergyCalculationStep = freeEnergyCalculationStep_ == step;
138 if (isEnergyCalculationStep || writeEnergy)
140 (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
141 [this, time, isEnergyCalculationStep, isFreeEnergyCalculationStep]() {
142 doStep(time, isEnergyCalculationStep, isFreeEnergyCalculationStep);
147 (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
148 [this]() { energyOutput_->recordNonEnergyStep(); }));
152 void EnergyElement::elementTeardown()
154 if (inputrec_->nstcalcenergy > 0 && isMasterRank_)
156 energyOutput_->printAverages(fplog_, groups_);
160 void EnergyElement::trajectoryWriterSetup(gmx_mdoutf* outf)
162 pull_t* pull_work = nullptr;
163 energyOutput_ = std::make_unique<EnergyOutput>(mdoutf_get_fp_ene(outf), top_global_, inputrec_,
164 pull_work, mdoutf_get_fp_dhdl(outf), false,
165 startingBehavior_, mdModulesNotifier_);
172 initializeEnergyHistory(startingBehavior_, observablesHistory_, energyOutput_.get());
174 // TODO: This probably doesn't really belong here...
175 // but we have all we need in this element,
176 // so we'll leave it here for now!
177 double io = compute_io(inputrec_, top_global_->natoms, *groups_, energyOutput_->numEnergyTerms(), 1);
178 if ((io > 2000) && isMasterRank_)
180 fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
182 if (!inputrec_->bContinuation)
184 real temp = enerd_->term[F_TEMP];
185 if (inputrec_->eI != eiVV)
187 /* Result of Ekin averaged over velocities of -half
188 * and +half step, while we only have -half step here.
192 fprintf(fplog_, "Initial temperature: %g K\n", temp);
196 ITrajectoryWriterCallbackPtr EnergyElement::registerTrajectoryWriterCallback(TrajectoryEvent event)
198 if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
200 return std::make_unique<ITrajectoryWriterCallback>(
201 [this](gmx_mdoutf* mdoutf, Step step, Time time, bool writeTrajectory,
202 bool writeLog) { write(mdoutf, step, time, writeTrajectory, writeLog); });
207 SignallerCallbackPtr EnergyElement::registerTrajectorySignallerCallback(gmx::TrajectoryEvent event)
209 if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
211 return std::make_unique<SignallerCallback>(
212 [this](Step step, Time /*unused*/) { energyWritingStep_ = step; });
217 SignallerCallbackPtr EnergyElement::registerEnergyCallback(EnergySignallerEvent event)
219 if (event == EnergySignallerEvent::EnergyCalculationStep && isMasterRank_)
221 return std::make_unique<SignallerCallback>(
222 [this](Step step, Time /*unused*/) { energyCalculationStep_ = step; });
224 if (event == EnergySignallerEvent::FreeEnergyCalculationStep && isMasterRank_)
226 return std::make_unique<SignallerCallback>(
227 [this](Step step, Time /*unused*/) { freeEnergyCalculationStep_ = step; });
232 void EnergyElement::doStep(Time time, bool isEnergyCalculationStep, bool isFreeEnergyCalculationStep)
234 enerd_->term[F_ETOT] = enerd_->term[F_EPOT] + enerd_->term[F_EKIN];
235 if (vRescaleThermostat_)
237 dummyLegacyState_.therm_integral = vRescaleThermostat_->thermostatIntegral();
239 if (freeEnergyPerturbationElement_)
241 sum_dhdl(enerd_, freeEnergyPerturbationElement_->constLambdaView(), *inputrec_->fepvals);
242 dummyLegacyState_.fep_state = freeEnergyPerturbationElement_->currentFEPState();
244 if (parrinelloRahmanBarostat_)
246 copy_mat(parrinelloRahmanBarostat_->boxVelocities(), dummyLegacyState_.boxv);
247 copy_mat(statePropagatorData_->constBox(), dummyLegacyState_.box);
249 if (integratorHasConservedEnergyQuantity(inputrec_))
251 enerd_->term[F_ECONSERVED] =
252 enerd_->term[F_ETOT] + NPT_energy(inputrec_, &dummyLegacyState_, nullptr);
254 energyOutput_->addDataAtEnergyStep(isFreeEnergyCalculationStep, isEnergyCalculationStep, time,
255 mdAtoms_->mdatoms()->tmass, enerd_, &dummyLegacyState_,
256 inputrec_->fepvals, inputrec_->expandedvals,
257 statePropagatorData_->constPreviousBox(), shakeVirial_,
258 forceVirial_, totalVirial_, pressure_, ekind_, muTot_, constr_);
261 void EnergyElement::write(gmx_mdoutf* outf, Step step, Time time, bool writeTrajectory, bool writeLog)
265 energyOutput_->printHeader(fplog_, step, time);
268 bool do_dr = do_per_step(step, inputrec_->nstdisreout);
269 bool do_or = do_per_step(step, inputrec_->nstorireout);
271 // energyOutput_->printAnnealingTemperatures(writeLog ? fplog_ : nullptr, groups_, &(inputrec_->opts));
273 energyOutput_->printStepToEnergyFile(mdoutf_get_fp_ene(outf), writeTrajectory, do_dr, do_or,
274 writeLog ? fplog_ : nullptr, step, time, fcd_, awh);
277 void EnergyElement::addToForceVirial(const tensor virial, Step step)
279 if (step > forceVirialStep_)
281 forceVirialStep_ = step;
282 clear_mat(forceVirial_);
284 m_add(forceVirial_, virial, forceVirial_);
287 void EnergyElement::addToConstraintVirial(const tensor virial, Step step)
289 if (step > shakeVirialStep_)
291 shakeVirialStep_ = step;
292 clear_mat(shakeVirial_);
294 m_add(shakeVirial_, virial, shakeVirial_);
297 rvec* EnergyElement::forceVirial(Step gmx_unused step)
299 if (step > forceVirialStep_)
301 forceVirialStep_ = step;
302 clear_mat(forceVirial_);
304 GMX_ASSERT(step >= forceVirialStep_ || forceVirialStep_ == -1,
305 "Asked for force virial of previous step.");
309 rvec* EnergyElement::constraintVirial(Step gmx_unused step)
311 if (step > shakeVirialStep_)
313 shakeVirialStep_ = step;
314 clear_mat(shakeVirial_);
316 GMX_ASSERT(step >= shakeVirialStep_ || shakeVirialStep_ == -1,
317 "Asked for constraint virial of previous step.");
321 rvec* EnergyElement::totalVirial(Step gmx_unused step)
323 if (step > totalVirialStep_)
325 totalVirialStep_ = step;
326 clear_mat(totalVirial_);
328 GMX_ASSERT(step >= totalVirialStep_ || totalVirialStep_ == -1,
329 "Asked for total virial of previous step.");
333 rvec* EnergyElement::pressure(Step gmx_unused step)
335 if (step > pressureStep_)
337 pressureStep_ = step;
338 clear_mat(pressure_);
340 GMX_ASSERT(step >= pressureStep_ || pressureStep_ == -1,
341 "Asked for pressure of previous step.");
345 real* EnergyElement::muTot()
350 gmx_enerdata_t* EnergyElement::enerdata()
355 gmx_ekindata_t* EnergyElement::ekindata()
360 bool* EnergyElement::needToSumEkinhOld()
362 return &needToSumEkinhOld_;
365 void EnergyElement::writeCheckpoint(t_state gmx_unused* localState, t_state* globalState)
369 if (needToSumEkinhOld_)
371 globalState->ekinstate.bUpToDate = false;
375 update_ekinstate(&globalState->ekinstate, ekind_);
376 globalState->ekinstate.bUpToDate = true;
378 energyOutput_->fillEnergyHistory(observablesHistory_->energyHistory.get());
382 void EnergyElement::initializeEnergyHistory(StartingBehavior startingBehavior,
383 ObservablesHistory* observablesHistory,
384 EnergyOutput* energyOutput)
386 if (startingBehavior != StartingBehavior::NewSimulation)
388 /* Restore from energy history if appending to output files */
389 if (startingBehavior == StartingBehavior::RestartWithAppending)
391 /* If no history is available (because a checkpoint is from before
392 * it was written) make a new one later, otherwise restore it.
394 if (observablesHistory->energyHistory)
396 energyOutput->restoreFromEnergyHistory(*observablesHistory->energyHistory);
399 else if (observablesHistory->energyHistory)
401 /* We might have read an energy history from checkpoint.
402 * As we are not appending, we want to restart the statistics.
403 * Free the allocated memory and reset the counts.
405 observablesHistory->energyHistory = {};
406 /* We might have read a pull history from checkpoint.
407 * We will still want to keep the statistics, so that the files
408 * can be joined and still be meaningful.
409 * This means that observablesHistory_->pullHistory
410 * should not be reset.
414 if (!observablesHistory->energyHistory)
416 observablesHistory->energyHistory = std::make_unique<energyhistory_t>();
418 if (!observablesHistory->pullHistory)
420 observablesHistory->pullHistory = std::make_unique<PullHistory>();
422 /* Set the initial energy history */
423 energyOutput->fillEnergyHistory(observablesHistory->energyHistory.get());
426 void EnergyElement::setVRescaleThermostat(const gmx::VRescaleThermostat* vRescaleThermostat)
428 vRescaleThermostat_ = vRescaleThermostat;
429 if (vRescaleThermostat_)
431 dummyLegacyState_.flags |= (1U << estTHERM_INT);
435 void EnergyElement::setParrinelloRahamnBarostat(const gmx::ParrinelloRahmanBarostat* parrinelloRahmanBarostat)
437 parrinelloRahmanBarostat_ = parrinelloRahmanBarostat;
438 if (parrinelloRahmanBarostat_)
440 dummyLegacyState_.flags |= (1U << estBOX) | (1U << estBOXV);