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