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