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