Use const ref for mtop and inputrec in EnergyOutput initialization
[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,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.
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>(mdoutf_get_fp_ene(outf),
171                                                    *top_global_,
172                                                    *inputrec_,
173                                                    pull_work,
174                                                    mdoutf_get_fp_dhdl(outf),
175                                                    false,
176                                                    startingBehavior_,
177                                                    simulationsShareState_,
178                                                    mdModulesNotifier_);
179
180     if (!isMasterRank_)
181     {
182         return;
183     }
184
185     initializeEnergyHistory(startingBehavior_, observablesHistory_, energyOutput_.get());
186
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_)
192     {
193         fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
194     }
195     if (!inputrec_->bContinuation)
196     {
197         real temp = enerd_->term[F_TEMP];
198         if (inputrec_->eI != IntegrationAlgorithm::VV)
199         {
200             /* Result of Ekin averaged over velocities of -half
201              * and +half step, while we only have -half step here.
202              */
203             temp *= 2;
204         }
205         fprintf(fplog_, "Initial temperature: %g K\n", temp);
206     }
207 }
208
209 std::optional<ITrajectoryWriterCallback> EnergyData::Element::registerTrajectoryWriterCallback(TrajectoryEvent event)
210 {
211     if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
212     {
213         return [this](gmx_mdoutf* mdoutf, Step step, Time time, bool writeTrajectory, bool writeLog) {
214             energyData_->write(mdoutf, step, time, writeTrajectory, writeLog);
215         };
216     }
217     return std::nullopt;
218 }
219
220 std::optional<SignallerCallback> EnergyData::Element::registerTrajectorySignallerCallback(gmx::TrajectoryEvent event)
221 {
222     if (event == TrajectoryEvent::EnergyWritingStep && isMasterRank_)
223     {
224         return [this](Step step, Time /*unused*/) { energyWritingStep_ = step; };
225     }
226     return std::nullopt;
227 }
228
229 std::optional<SignallerCallback> EnergyData::Element::registerEnergyCallback(EnergySignallerEvent event)
230 {
231     if (event == EnergySignallerEvent::EnergyCalculationStep && isMasterRank_)
232     {
233         return [this](Step step, Time /*unused*/) { energyCalculationStep_ = step; };
234     }
235     if (event == EnergySignallerEvent::FreeEnergyCalculationStep && isMasterRank_)
236     {
237         return [this](Step step, Time /*unused*/) { freeEnergyCalculationStep_ = step; };
238     }
239     return std::nullopt;
240 }
241
242 void EnergyData::doStep(Time time, bool isEnergyCalculationStep, bool isFreeEnergyCalculationStep)
243 {
244     enerd_->term[F_ETOT] = enerd_->term[F_EPOT] + enerd_->term[F_EKIN];
245     if (freeEnergyPerturbationData_)
246     {
247         accumulateKineticLambdaComponents(
248                 enerd_, freeEnergyPerturbationData_->constLambdaView(), *inputrec_->fepvals);
249     }
250     if (integratorHasConservedEnergyQuantity(inputrec_))
251     {
252         enerd_->term[F_ECONSERVED] =
253                 enerd_->term[F_ETOT]
254                 + (velocityScalingTemperatureCoupling_
255                            ? velocityScalingTemperatureCoupling_->conservedEnergyContribution()
256                            : 0)
257                 + (parrinelloRahmanBarostat_ ? parrinelloRahmanBarostat_->conservedEnergyContribution() : 0);
258     }
259     matrix nullMatrix = {};
260     energyOutput_->addDataAtEnergyStep(
261             isFreeEnergyCalculationStep,
262             isEnergyCalculationStep,
263             time,
264             mdAtoms_->mdatoms()->tmass,
265             enerd_,
266             inputrec_->fepvals.get(),
267             inputrec_->expandedvals.get(),
268             statePropagatorData_->constPreviousBox(),
269             PTCouplingArrays({ parrinelloRahmanBarostat_ ? parrinelloRahmanBarostat_->boxVelocities() : nullMatrix,
270                                {},
271                                {},
272                                {},
273                                {} }),
274             freeEnergyPerturbationData_ ? freeEnergyPerturbationData_->currentFEPState() : 0,
275             shakeVirial_,
276             forceVirial_,
277             totalVirial_,
278             pressure_,
279             ekind_,
280             muTot_,
281             constr_);
282 }
283
284 void EnergyData::write(gmx_mdoutf* outf, Step step, Time time, bool writeTrajectory, bool writeLog)
285 {
286     if (writeLog)
287     {
288         energyOutput_->printHeader(fplog_, step, time);
289     }
290
291     bool do_dr = do_per_step(step, inputrec_->nstdisreout);
292     bool do_or = do_per_step(step, inputrec_->nstorireout);
293
294     // energyOutput_->printAnnealingTemperatures(writeLog ? fplog_ : nullptr, groups_, &(inputrec_->opts));
295     Awh* awh = nullptr;
296     energyOutput_->printStepToEnergyFile(
297             mdoutf_get_fp_ene(outf), writeTrajectory, do_dr, do_or, writeLog ? fplog_ : nullptr, step, time, fcd_, awh);
298 }
299
300 void EnergyData::addToForceVirial(const tensor virial, Step step)
301 {
302     if (step > forceVirialStep_)
303     {
304         forceVirialStep_ = step;
305         clear_mat(forceVirial_);
306     }
307     m_add(forceVirial_, virial, forceVirial_);
308 }
309
310 void EnergyData::addToConstraintVirial(const tensor virial, Step step)
311 {
312     if (step > shakeVirialStep_)
313     {
314         shakeVirialStep_ = step;
315         clear_mat(shakeVirial_);
316     }
317     m_add(shakeVirial_, virial, shakeVirial_);
318 }
319
320 rvec* EnergyData::forceVirial(Step gmx_unused step)
321 {
322     if (step > forceVirialStep_)
323     {
324         forceVirialStep_ = step;
325         clear_mat(forceVirial_);
326     }
327     GMX_ASSERT(step >= forceVirialStep_ || forceVirialStep_ == -1,
328                "Asked for force virial of previous step.");
329     return forceVirial_;
330 }
331
332 rvec* EnergyData::constraintVirial(Step gmx_unused step)
333 {
334     if (step > shakeVirialStep_)
335     {
336         shakeVirialStep_ = step;
337         clear_mat(shakeVirial_);
338     }
339     GMX_ASSERT(step >= shakeVirialStep_ || shakeVirialStep_ == -1,
340                "Asked for constraint virial of previous step.");
341     return shakeVirial_;
342 }
343
344 rvec* EnergyData::totalVirial(Step gmx_unused step)
345 {
346     if (step > totalVirialStep_)
347     {
348         totalVirialStep_ = step;
349         clear_mat(totalVirial_);
350     }
351     GMX_ASSERT(step >= totalVirialStep_ || totalVirialStep_ == -1,
352                "Asked for total virial of previous step.");
353     return totalVirial_;
354 }
355
356 rvec* EnergyData::pressure(Step gmx_unused step)
357 {
358     if (step > pressureStep_)
359     {
360         pressureStep_ = step;
361         clear_mat(pressure_);
362     }
363     GMX_ASSERT(step >= pressureStep_ || pressureStep_ == -1, "Asked for pressure of previous step.");
364     return pressure_;
365 }
366
367 real* EnergyData::muTot()
368 {
369     return muTot_;
370 }
371
372 gmx_enerdata_t* EnergyData::enerdata()
373 {
374     return enerd_;
375 }
376
377 gmx_ekindata_t* EnergyData::ekindata()
378 {
379     return ekind_;
380 }
381
382 bool* EnergyData::needToSumEkinhOld()
383 {
384     return &needToSumEkinhOld_;
385 }
386
387 bool EnergyData::hasReadEkinFromCheckpoint() const
388 {
389     return hasReadEkinFromCheckpoint_;
390 }
391
392 namespace
393 {
394 /*!
395  * \brief Enum describing the contents EnergyData::Element writes to modular checkpoint
396  *
397  * When changing the checkpoint content, add a new element just above Count, and adjust the
398  * checkpoint functionality.
399  */
400 enum class CheckpointVersion
401 {
402     Base, //!< First version of modular checkpointing
403     Count //!< Number of entries. Add new versions right above this!
404 };
405 constexpr auto c_currentVersion = CheckpointVersion(int(CheckpointVersion::Count) - 1);
406 } // namespace
407
408 template<CheckpointDataOperation operation>
409 void EnergyData::Element::doCheckpointData(CheckpointData<operation>* checkpointData)
410 {
411     checkpointVersion(checkpointData, "EnergyData version", c_currentVersion);
412
413     energyData_->observablesHistory_->energyHistory->doCheckpoint<operation>(
414             checkpointData->subCheckpointData("energy history"));
415     energyData_->ekinstate_.doCheckpoint<operation>(checkpointData->subCheckpointData("ekinstate"));
416 }
417
418 void EnergyData::Element::saveCheckpointState(std::optional<WriteCheckpointData> checkpointData,
419                                               const t_commrec*                   cr)
420 {
421     if (MASTER(cr))
422     {
423         if (energyData_->needToSumEkinhOld_)
424         {
425             energyData_->ekinstate_.bUpToDate = false;
426         }
427         else
428         {
429             update_ekinstate(&energyData_->ekinstate_, energyData_->ekind_);
430             energyData_->ekinstate_.bUpToDate = true;
431         }
432         energyData_->energyOutput_->fillEnergyHistory(
433                 energyData_->observablesHistory_->energyHistory.get());
434         doCheckpointData<CheckpointDataOperation::Write>(&checkpointData.value());
435     }
436 }
437
438 void EnergyData::Element::restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData,
439                                                  const t_commrec*                  cr)
440 {
441     if (MASTER(cr))
442     {
443         doCheckpointData<CheckpointDataOperation::Read>(&checkpointData.value());
444     }
445     energyData_->hasReadEkinFromCheckpoint_ = MASTER(cr) ? energyData_->ekinstate_.bUpToDate : false;
446     if (PAR(cr))
447     {
448         gmx_bcast(sizeof(hasReadEkinFromCheckpoint_),
449                   &energyData_->hasReadEkinFromCheckpoint_,
450                   cr->mpi_comm_mygroup);
451     }
452     if (energyData_->hasReadEkinFromCheckpoint_)
453     {
454         // this takes care of broadcasting from master to agents
455         restore_ekinstate_from_state(cr, energyData_->ekind_, &energyData_->ekinstate_);
456     }
457 }
458
459 const std::string& EnergyData::Element::clientID()
460 {
461     return identifier_;
462 }
463
464 void EnergyData::initializeEnergyHistory(StartingBehavior    startingBehavior,
465                                          ObservablesHistory* observablesHistory,
466                                          EnergyOutput*       energyOutput)
467 {
468     if (startingBehavior != StartingBehavior::NewSimulation)
469     {
470         /* Restore from energy history if appending to output files */
471         if (startingBehavior == StartingBehavior::RestartWithAppending)
472         {
473             /* If no history is available (because a checkpoint is from before
474              * it was written) make a new one later, otherwise restore it.
475              */
476             if (observablesHistory->energyHistory)
477             {
478                 energyOutput->restoreFromEnergyHistory(*observablesHistory->energyHistory);
479             }
480         }
481         else if (observablesHistory->energyHistory)
482         {
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.
486              */
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.
493              */
494         }
495     }
496     if (!observablesHistory->energyHistory)
497     {
498         observablesHistory->energyHistory = std::make_unique<energyhistory_t>();
499     }
500     if (!observablesHistory->pullHistory)
501     {
502         observablesHistory->pullHistory = std::make_unique<PullHistory>();
503     }
504     /* Set the initial energy history */
505     energyOutput->fillEnergyHistory(observablesHistory->energyHistory.get());
506 }
507
508 void EnergyData::setVelocityScalingTemperatureCoupling(const VelocityScalingTemperatureCoupling* velocityScalingTemperatureCoupling)
509 {
510     velocityScalingTemperatureCoupling_ = velocityScalingTemperatureCoupling;
511 }
512
513 void EnergyData::setParrinelloRahamnBarostat(const gmx::ParrinelloRahmanBarostat* parrinelloRahmanBarostat)
514 {
515     parrinelloRahmanBarostat_ = parrinelloRahmanBarostat;
516 }
517
518 EnergyData::Element* EnergyData::element()
519 {
520     return element_.get();
521 }
522
523 EnergyData::Element::Element(EnergyData* energyData, bool isMasterRank) :
524     energyData_(energyData),
525     isMasterRank_(isMasterRank),
526     energyWritingStep_(-1),
527     energyCalculationStep_(-1),
528     freeEnergyCalculationStep_(-1)
529 {
530 }
531
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)
539 {
540     return energyData->element();
541 }
542
543 } // namespace gmx