Print energy conservation in modular simulator
[alexxy/gromacs.git] / src / gromacs / modularsimulator / energydata.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief Defines the microstate for the modular simulator
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_modularsimulator
40  */
41
42 #include "gmxpre.h"
43
44 #include "energydata.h"
45
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"
66
67 #include "freeenergyperturbationdata.h"
68 #include "modularsimulator.h"
69 #include "parrinellorahmanbarostat.h"
70 #include "simulatoralgorithm.h"
71 #include "statepropagatordata.h"
72 #include "velocityscalingtemperaturecoupling.h"
73
74 struct pull_t;
75 class t_state;
76
77 namespace gmx
78 {
79 class Awh;
80
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,
89                        FILE*                       fplog,
90                        t_fcdata*                   fcd,
91                        const MdModulesNotifier&    mdModulesNotifier,
92                        bool                        isMasterRank,
93                        ObservablesHistory*         observablesHistory,
94                        StartingBehavior            startingBehavior,
95                        bool                        simulationsShareState) :
96     element_(std::make_unique<Element>(this, isMasterRank)),
97     isMasterRank_(isMasterRank),
98     forceVirialStep_(-1),
99     shakeVirialStep_(-1),
100     totalVirialStep_(-1),
101     pressureStep_(-1),
102     needToSumEkinhOld_(false),
103     hasReadEkinFromCheckpoint_(false),
104     startingBehavior_(startingBehavior),
105     statePropagatorData_(statePropagatorData),
106     freeEnergyPerturbationData_(freeEnergyPerturbationData),
107     velocityScalingTemperatureCoupling_(nullptr),
108     parrinelloRahmanBarostat_(nullptr),
109     inputrec_(inputrec),
110     top_global_(globalTopology),
111     mdAtoms_(mdAtoms),
112     enerd_(enerd),
113     ekind_(ekind),
114     constr_(constr),
115     fplog_(fplog),
116     fcd_(fcd),
117     mdModulesNotifier_(mdModulesNotifier),
118     groups_(&globalTopology->groups),
119     observablesHistory_(observablesHistory),
120     simulationsShareState_(simulationsShareState)
121 {
122     clear_mat(forceVirial_);
123     clear_mat(shakeVirial_);
124     clear_mat(totalVirial_);
125     clear_mat(pressure_);
126     clear_rvec(muTot_);
127
128     init_ekinstate(&ekinstate_, inputrec_);
129     observablesHistory_->energyHistory = std::make_unique<energyhistory_t>();
130 }
131
132 void EnergyData::Element::scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction)
133 {
134     if (!isMasterRank_)
135     {
136         return;
137     }
138     auto writeEnergy                 = energyWritingStep_ == step;
139     auto isEnergyCalculationStep     = energyCalculationStep_ == step;
140     auto isFreeEnergyCalculationStep = freeEnergyCalculationStep_ == step;
141     if (isEnergyCalculationStep || writeEnergy)
142     {
143         registerRunFunction([this, time, isEnergyCalculationStep, isFreeEnergyCalculationStep]() {
144             energyData_->doStep(time, isEnergyCalculationStep, isFreeEnergyCalculationStep);
145         });
146     }
147     else
148     {
149         registerRunFunction([this]() { energyData_->energyOutput_->recordNonEnergyStep(); });
150     }
151 }
152
153 void EnergyData::teardown()
154 {
155     if (inputrec_->nstcalcenergy > 0 && isMasterRank_)
156     {
157         energyOutput_->printEnergyConservation(fplog_, inputrec_->simulation_part, EI_MD(inputrec_->eI));
158         energyOutput_->printAverages(fplog_, groups_);
159     }
160 }
161
162 void EnergyData::Element::trajectoryWriterSetup(gmx_mdoutf* outf)
163 {
164     energyData_->setup(outf);
165 }
166
167 void EnergyData::setup(gmx_mdoutf* outf)
168 {
169     pull_t* pull_work = nullptr;
170     energyOutput_     = std::make_unique<EnergyOutput>(
171             mdoutf_get_fp_ene(outf), top_global_, inputrec_, pull_work, mdoutf_get_fp_dhdl(outf),
172             false, startingBehavior_, simulationsShareState_, mdModulesNotifier_);
173
174     if (!isMasterRank_)
175     {
176         return;
177     }
178
179     initializeEnergyHistory(startingBehavior_, observablesHistory_, energyOutput_.get());
180
181     // TODO: This probably doesn't really belong here...
182     //       but we have all we need in this element,
183     //       so we'll leave it here for now!
184     double io = compute_io(inputrec_, top_global_->natoms, *groups_, energyOutput_->numEnergyTerms(), 1);
185     if ((io > 2000) && isMasterRank_)
186     {
187         fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
188     }
189     if (!inputrec_->bContinuation)
190     {
191         real temp = enerd_->term[F_TEMP];
192         if (inputrec_->eI != eiVV)
193         {
194             /* Result of Ekin averaged over velocities of -half
195              * and +half step, while we only have -half step here.
196              */
197             temp *= 2;
198         }
199         fprintf(fplog_, "Initial temperature: %g K\n", temp);
200     }
201 }
202
203 std::optional<ITrajectoryWriterCallback> EnergyData::Element::registerTrajectoryWriterCallback(TrajectoryEvent event)
204 {
205     if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
206     {
207         return [this](gmx_mdoutf* mdoutf, Step step, Time time, bool writeTrajectory, bool writeLog) {
208             energyData_->write(mdoutf, step, time, writeTrajectory, writeLog);
209         };
210     }
211     return std::nullopt;
212 }
213
214 std::optional<SignallerCallback> EnergyData::Element::registerTrajectorySignallerCallback(gmx::TrajectoryEvent event)
215 {
216     if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
217     {
218         return [this](Step step, Time /*unused*/) { energyWritingStep_ = step; };
219     }
220     return std::nullopt;
221 }
222
223 std::optional<SignallerCallback> EnergyData::Element::registerEnergyCallback(EnergySignallerEvent event)
224 {
225     if (event == EnergySignallerEvent::EnergyCalculationStep && isMasterRank_)
226     {
227         return [this](Step step, Time /*unused*/) { energyCalculationStep_ = step; };
228     }
229     if (event == EnergySignallerEvent::FreeEnergyCalculationStep && isMasterRank_)
230     {
231         return [this](Step step, Time /*unused*/) { freeEnergyCalculationStep_ = step; };
232     }
233     return std::nullopt;
234 }
235
236 void EnergyData::doStep(Time time, bool isEnergyCalculationStep, bool isFreeEnergyCalculationStep)
237 {
238     enerd_->term[F_ETOT] = enerd_->term[F_EPOT] + enerd_->term[F_EKIN];
239     if (freeEnergyPerturbationData_)
240     {
241         accumulateKineticLambdaComponents(enerd_, freeEnergyPerturbationData_->constLambdaView(),
242                                           *inputrec_->fepvals);
243     }
244     if (integratorHasConservedEnergyQuantity(inputrec_))
245     {
246         enerd_->term[F_ECONSERVED] =
247                 enerd_->term[F_ETOT]
248                 + (velocityScalingTemperatureCoupling_
249                            ? velocityScalingTemperatureCoupling_->conservedEnergyContribution()
250                            : 0)
251                 + (parrinelloRahmanBarostat_ ? parrinelloRahmanBarostat_->conservedEnergyContribution() : 0);
252     }
253     matrix nullMatrix = {};
254     energyOutput_->addDataAtEnergyStep(
255             isFreeEnergyCalculationStep, isEnergyCalculationStep, time, mdAtoms_->mdatoms()->tmass, enerd_,
256             inputrec_->fepvals, inputrec_->expandedvals, statePropagatorData_->constPreviousBox(),
257             PTCouplingArrays({ parrinelloRahmanBarostat_ ? parrinelloRahmanBarostat_->boxVelocities() : nullMatrix,
258                                {},
259                                {},
260                                {},
261                                {} }),
262             freeEnergyPerturbationData_ ? freeEnergyPerturbationData_->currentFEPState() : 0,
263             shakeVirial_, forceVirial_, totalVirial_, pressure_, ekind_, muTot_, constr_);
264 }
265
266 void EnergyData::write(gmx_mdoutf* outf, Step step, Time time, bool writeTrajectory, bool writeLog)
267 {
268     if (writeLog)
269     {
270         energyOutput_->printHeader(fplog_, step, time);
271     }
272
273     bool do_dr = do_per_step(step, inputrec_->nstdisreout);
274     bool do_or = do_per_step(step, inputrec_->nstorireout);
275
276     // energyOutput_->printAnnealingTemperatures(writeLog ? fplog_ : nullptr, groups_, &(inputrec_->opts));
277     Awh* awh = nullptr;
278     energyOutput_->printStepToEnergyFile(mdoutf_get_fp_ene(outf), writeTrajectory, do_dr, do_or,
279                                          writeLog ? fplog_ : nullptr, step, time, fcd_, awh);
280 }
281
282 void EnergyData::addToForceVirial(const tensor virial, Step step)
283 {
284     if (step > forceVirialStep_)
285     {
286         forceVirialStep_ = step;
287         clear_mat(forceVirial_);
288     }
289     m_add(forceVirial_, virial, forceVirial_);
290 }
291
292 void EnergyData::addToConstraintVirial(const tensor virial, Step step)
293 {
294     if (step > shakeVirialStep_)
295     {
296         shakeVirialStep_ = step;
297         clear_mat(shakeVirial_);
298     }
299     m_add(shakeVirial_, virial, shakeVirial_);
300 }
301
302 rvec* EnergyData::forceVirial(Step gmx_unused step)
303 {
304     if (step > forceVirialStep_)
305     {
306         forceVirialStep_ = step;
307         clear_mat(forceVirial_);
308     }
309     GMX_ASSERT(step >= forceVirialStep_ || forceVirialStep_ == -1,
310                "Asked for force virial of previous step.");
311     return forceVirial_;
312 }
313
314 rvec* EnergyData::constraintVirial(Step gmx_unused step)
315 {
316     if (step > shakeVirialStep_)
317     {
318         shakeVirialStep_ = step;
319         clear_mat(shakeVirial_);
320     }
321     GMX_ASSERT(step >= shakeVirialStep_ || shakeVirialStep_ == -1,
322                "Asked for constraint virial of previous step.");
323     return shakeVirial_;
324 }
325
326 rvec* EnergyData::totalVirial(Step gmx_unused step)
327 {
328     if (step > totalVirialStep_)
329     {
330         totalVirialStep_ = step;
331         clear_mat(totalVirial_);
332     }
333     GMX_ASSERT(step >= totalVirialStep_ || totalVirialStep_ == -1,
334                "Asked for total virial of previous step.");
335     return totalVirial_;
336 }
337
338 rvec* EnergyData::pressure(Step gmx_unused step)
339 {
340     if (step > pressureStep_)
341     {
342         pressureStep_ = step;
343         clear_mat(pressure_);
344     }
345     GMX_ASSERT(step >= pressureStep_ || pressureStep_ == -1,
346                "Asked for pressure of previous step.");
347     return pressure_;
348 }
349
350 real* EnergyData::muTot()
351 {
352     return muTot_;
353 }
354
355 gmx_enerdata_t* EnergyData::enerdata()
356 {
357     return enerd_;
358 }
359
360 gmx_ekindata_t* EnergyData::ekindata()
361 {
362     return ekind_;
363 }
364
365 bool* EnergyData::needToSumEkinhOld()
366 {
367     return &needToSumEkinhOld_;
368 }
369
370 bool EnergyData::hasReadEkinFromCheckpoint() const
371 {
372     return hasReadEkinFromCheckpoint_;
373 }
374
375 namespace
376 {
377 /*!
378  * \brief Enum describing the contents EnergyData::Element writes to modular checkpoint
379  *
380  * When changing the checkpoint content, add a new element just above Count, and adjust the
381  * checkpoint functionality.
382  */
383 enum class CheckpointVersion
384 {
385     Base, //!< First version of modular checkpointing
386     Count //!< Number of entries. Add new versions right above this!
387 };
388 constexpr auto c_currentVersion = CheckpointVersion(int(CheckpointVersion::Count) - 1);
389 } // namespace
390
391 template<CheckpointDataOperation operation>
392 void EnergyData::Element::doCheckpointData(CheckpointData<operation>* checkpointData)
393 {
394     checkpointVersion(checkpointData, "EnergyData version", c_currentVersion);
395
396     energyData_->observablesHistory_->energyHistory->doCheckpoint<operation>(
397             checkpointData->subCheckpointData("energy history"));
398     energyData_->ekinstate_.doCheckpoint<operation>(checkpointData->subCheckpointData("ekinstate"));
399 }
400
401 void EnergyData::Element::saveCheckpointState(std::optional<WriteCheckpointData> checkpointData,
402                                               const t_commrec*                   cr)
403 {
404     if (MASTER(cr))
405     {
406         if (energyData_->needToSumEkinhOld_)
407         {
408             energyData_->ekinstate_.bUpToDate = false;
409         }
410         else
411         {
412             update_ekinstate(&energyData_->ekinstate_, energyData_->ekind_);
413             energyData_->ekinstate_.bUpToDate = true;
414         }
415         energyData_->energyOutput_->fillEnergyHistory(
416                 energyData_->observablesHistory_->energyHistory.get());
417         doCheckpointData<CheckpointDataOperation::Write>(&checkpointData.value());
418     }
419 }
420
421 void EnergyData::Element::restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData,
422                                                  const t_commrec*                  cr)
423 {
424     if (MASTER(cr))
425     {
426         doCheckpointData<CheckpointDataOperation::Read>(&checkpointData.value());
427     }
428     energyData_->hasReadEkinFromCheckpoint_ = MASTER(cr) ? energyData_->ekinstate_.bUpToDate : false;
429     if (PAR(cr))
430     {
431         gmx_bcast(sizeof(hasReadEkinFromCheckpoint_), &energyData_->hasReadEkinFromCheckpoint_,
432                   cr->mpi_comm_mygroup);
433     }
434     if (energyData_->hasReadEkinFromCheckpoint_)
435     {
436         // this takes care of broadcasting from master to agents
437         restore_ekinstate_from_state(cr, energyData_->ekind_, &energyData_->ekinstate_);
438     }
439 }
440
441 const std::string& EnergyData::Element::clientID()
442 {
443     return identifier_;
444 }
445
446 void EnergyData::initializeEnergyHistory(StartingBehavior    startingBehavior,
447                                          ObservablesHistory* observablesHistory,
448                                          EnergyOutput*       energyOutput)
449 {
450     if (startingBehavior != StartingBehavior::NewSimulation)
451     {
452         /* Restore from energy history if appending to output files */
453         if (startingBehavior == StartingBehavior::RestartWithAppending)
454         {
455             /* If no history is available (because a checkpoint is from before
456              * it was written) make a new one later, otherwise restore it.
457              */
458             if (observablesHistory->energyHistory)
459             {
460                 energyOutput->restoreFromEnergyHistory(*observablesHistory->energyHistory);
461             }
462         }
463         else if (observablesHistory->energyHistory)
464         {
465             /* We might have read an energy history from checkpoint.
466              * As we are not appending, we want to restart the statistics.
467              * Free the allocated memory and reset the counts.
468              */
469             observablesHistory->energyHistory = {};
470             /* We might have read a pull history from checkpoint.
471              * We will still want to keep the statistics, so that the files
472              * can be joined and still be meaningful.
473              * This means that observablesHistory_->pullHistory
474              * should not be reset.
475              */
476         }
477     }
478     if (!observablesHistory->energyHistory)
479     {
480         observablesHistory->energyHistory = std::make_unique<energyhistory_t>();
481     }
482     if (!observablesHistory->pullHistory)
483     {
484         observablesHistory->pullHistory = std::make_unique<PullHistory>();
485     }
486     /* Set the initial energy history */
487     energyOutput->fillEnergyHistory(observablesHistory->energyHistory.get());
488 }
489
490 void EnergyData::setVelocityScalingTemperatureCoupling(const VelocityScalingTemperatureCoupling* velocityScalingTemperatureCoupling)
491 {
492     velocityScalingTemperatureCoupling_ = velocityScalingTemperatureCoupling;
493 }
494
495 void EnergyData::setParrinelloRahamnBarostat(const gmx::ParrinelloRahmanBarostat* parrinelloRahmanBarostat)
496 {
497     parrinelloRahmanBarostat_ = parrinelloRahmanBarostat;
498 }
499
500 EnergyData::Element* EnergyData::element()
501 {
502     return element_.get();
503 }
504
505 EnergyData::Element::Element(EnergyData* energyData, bool isMasterRank) :
506     energyData_(energyData),
507     isMasterRank_(isMasterRank),
508     energyWritingStep_(-1),
509     energyCalculationStep_(-1),
510     freeEnergyCalculationStep_(-1)
511 {
512 }
513
514 ISimulatorElement* EnergyData::Element::getElementPointerImpl(
515         LegacySimulatorData gmx_unused*        legacySimulatorData,
516         ModularSimulatorAlgorithmBuilderHelper gmx_unused* builderHelper,
517         StatePropagatorData gmx_unused* statePropagatorData,
518         EnergyData*                     energyData,
519         FreeEnergyPerturbationData gmx_unused* freeEnergyPerturbationData,
520         GlobalCommunicationHelper gmx_unused* globalCommunicationHelper)
521 {
522     return energyData->element();
523 }
524
525 } // namespace gmx