Introduce plumbing for ObservablesReducer
[alexxy/gromacs.git] / src / gromacs / modularsimulator / energydata.h
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 Declares the energy element for the modular simulator
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_modularsimulator
40  *
41  * This header is only used within the modular simulator module
42  */
43
44 #ifndef GMX_ENERGYELEMENT_MICROSTATE_H
45 #define GMX_ENERGYELEMENT_MICROSTATE_H
46
47 #include "gromacs/mdtypes/state.h"
48
49 #include "modularsimulatorinterfaces.h"
50
51 class gmx_ekindata_t;
52 struct gmx_enerdata_t;
53 struct gmx_mtop_t;
54 struct ObservablesHistory;
55 struct t_fcdata;
56 struct t_inputrec;
57 struct SimulationGroups;
58
59 namespace gmx
60 {
61 enum class StartingBehavior;
62 class Constraints;
63 class EnergyOutput;
64 class FreeEnergyPerturbationData;
65 class GlobalCommunicationHelper;
66 class LegacySimulatorData;
67 class MDAtoms;
68 class ModularSimulatorAlgorithmBuilderHelper;
69 class ObservablesReducer;
70 class ParrinelloRahmanBarostat;
71 class StatePropagatorData;
72 class VelocityScalingTemperatureCoupling;
73 struct MDModulesNotifiers;
74
75 //! Function type for elements contributing energy
76 using EnergyContribution = std::function<real(Step, Time)>;
77
78 /*! \internal
79  * \ingroup module_modularsimulator
80  * \brief Data class managing energies
81  *
82  * The EnergyData owns the EnergyObject,
83  * the tensors for the different virials and the pressure as well as
84  * the total dipole vector. It has a member class which is part of the
85  * simulator loop and and is responsible
86  * for saving energy data and writing it to trajectory.
87  *
88  * The EnergyData offers an interface to add virial contributions,
89  * but also allows access to the raw pointers to tensor data, the
90  * dipole vector, and the legacy energy data structures.
91  *
92  * The EnergyData owns an object of type EnergyData::Element,
93  * which takes part in the simulation loop, allowing to record
94  * and output energies during the simulation.
95  */
96 class EnergyData final
97 {
98 public:
99     //! Constructor
100     EnergyData(StatePropagatorData*        statePropagatorData,
101                FreeEnergyPerturbationData* freeEnergyPerturbationData,
102                const gmx_mtop_t&           globalTopology,
103                const t_inputrec*           inputrec,
104                const MDAtoms*              mdAtoms,
105                gmx_enerdata_t*             enerd,
106                gmx_ekindata_t*             ekind,
107                const Constraints*          constr,
108                FILE*                       fplog,
109                t_fcdata*                   fcd,
110                const MDModulesNotifiers&   mdModulesNotifiers,
111                bool                        isMasterRank,
112                ObservablesHistory*         observablesHistory,
113                StartingBehavior            startingBehavior,
114                bool                        simulationsShareState);
115
116     /*! \brief Final output
117      *
118      * Prints the averages to log. This is called from ModularSimulatorAlgorithm.
119      *
120      * \see ModularSimulatorAlgorithm::teardown
121      */
122     void teardown();
123
124     /*! \brief Add contribution to force virial
125      *
126      * This automatically resets the tensor if the step is higher
127      * than the current step, starting the tensor calculation for
128      * a new step at zero. Otherwise, it adds the new contribution
129      * to the existing virial.
130      */
131     void addToForceVirial(const tensor virial, Step step);
132
133     /*! \brief Add contribution to constraint virial
134      *
135      * This automatically resets the tensor if the step is higher
136      * than the current step, starting the tensor calculation for
137      * a new step at zero. Otherwise, it adds the new contribution
138      * to the existing virial.
139      */
140     void addToConstraintVirial(const tensor virial, Step step);
141
142     /*! \brief Get pointer to force virial tensor
143      *
144      * Allows access to the raw pointer to the tensor.
145      */
146     rvec* forceVirial(Step step);
147
148     /*! \brief Get pointer to constraint virial tensor
149      *
150      * Allows access to the raw pointer to the tensor.
151      */
152     rvec* constraintVirial(Step step);
153
154     /*! \brief Get pointer to total virial tensor
155      *
156      * Allows access to the raw pointer to the tensor.
157      */
158     rvec* totalVirial(Step step);
159
160     /*! \brief Get pointer to pressure tensor
161      *
162      * Allows access to the raw pointer to the tensor.
163      */
164     rvec* pressure(Step step);
165
166     /*! \brief Get pointer to mu_tot
167      *
168      * Allows access to the raw pointer to the dipole vector.
169      */
170     real* muTot();
171
172     /*! \brief Get pointer to energy structure
173      *
174      */
175     gmx_enerdata_t* enerdata();
176
177     /*! \brief Get const pointer to energy structure
178      *
179      */
180     const gmx_enerdata_t* enerdata() const;
181
182     /*! \brief Get pointer to kinetic energy structure
183      *
184      */
185     gmx_ekindata_t* ekindata();
186
187     /*! \brief Get pointer to needToSumEkinhOld
188      *
189      */
190     bool* needToSumEkinhOld();
191
192     /*! \brief Whether kinetic energy was read from checkpoint
193      *
194      * This is needed by the compute globals element
195      * TODO: Remove this when moving global reduction to client system (#3421)
196      */
197     [[nodiscard]] bool hasReadEkinFromCheckpoint() const;
198
199     /*! \brief Add conserved energy contribution
200      *
201      * This allows other elements to register callbacks for contributions to
202      * the conserved energy term.
203      */
204     void addConservedEnergyContribution(EnergyContribution&& energyContribution);
205
206     /*! \brief set Parrinello-Rahman barostat
207      *
208      * This allows to set a pointer to the Parrinello-Rahman barostat used to
209      * print the box velocities.
210      */
211     void setParrinelloRahmanBoxVelocities(std::function<const rvec*()>&& parrinelloRahmanBoxVelocities);
212
213     /*! \brief Initialize energy history
214      *
215      * Kept as a static function to allow usage from legacy code
216      * \todo Make member function once legacy use is not needed anymore
217      */
218     static void initializeEnergyHistory(StartingBehavior    startingBehavior,
219                                         ObservablesHistory* observablesHistory,
220                                         EnergyOutput*       energyOutput);
221
222     /*! \brief Request (local) kinetic energy update
223      */
224     void updateKineticEnergy();
225
226     //! The element taking part in the simulator loop
227     class Element;
228     //! Get pointer to element (whose lifetime is managed by this)
229     Element* element();
230
231 private:
232     /*! \brief Setup (needs file pointer)
233      *
234      * Initializes the EnergyOutput object, and does some logging output.
235      *
236      * \param mdoutf  File pointer
237      */
238     void setup(gmx_mdoutf* mdoutf);
239
240     /*! \brief Save data at energy steps
241      *
242      * \param step  The current step
243      * \param time  The current time
244      * \param isEnergyCalculationStep  Whether the current step is an energy calculation step
245      * \param isFreeEnergyCalculationStep  Whether the current step is a free energy calculation step
246      */
247     void doStep(Step step, Time time, bool isEnergyCalculationStep, bool isFreeEnergyCalculationStep);
248
249     /*! \brief Write to energy trajectory
250      *
251      * This is only called by master - writes energy to trajectory and to log.
252      */
253     void write(gmx_mdoutf* outf, Step step, Time time, bool writeTrajectory, bool writeLog);
254
255     /*
256      * Data owned by EnergyData
257      */
258     //! The element
259     std::unique_ptr<Element> element_;
260     //! The energy output object
261     std::unique_ptr<EnergyOutput> energyOutput_;
262     //! Helper object to checkpoint kinetic energy data
263     ekinstate_t ekinstate_;
264
265     //! Whether this is the master rank
266     const bool isMasterRank_;
267
268     //! The force virial tensor
269     tensor forceVirial_;
270     //! The constraint virial tensor
271     tensor shakeVirial_;
272     //! The total virial tensor
273     tensor totalVirial_;
274     //! The pressure tensor
275     tensor pressure_;
276     //! The total dipole moment
277     rvec muTot_;
278
279     //! The step number of the current force virial tensor
280     Step forceVirialStep_;
281     //! The step number of the current constraint virial tensor
282     Step shakeVirialStep_;
283     //! The step number of the current total virial tensor
284     Step totalVirialStep_;
285     //! The step number of the current pressure tensor
286     Step pressureStep_;
287
288     //! Whether ekinh_old needs to be summed up (set by compute globals)
289     bool needToSumEkinhOld_;
290     //! Whether we have read ekin from checkpoint
291     bool hasReadEkinFromCheckpoint_;
292
293     //! Describes how the simulation (re)starts
294     const StartingBehavior startingBehavior_;
295
296     /*
297      * Pointers to Simulator data
298      */
299     // TODO: Clarify relationship to data objects and find a more robust alternative to raw pointers (#3583)
300     //! Pointer to the state propagator data
301     StatePropagatorData* statePropagatorData_;
302     //! Pointer to the free energy perturbation data
303     FreeEnergyPerturbationData* freeEnergyPerturbationData_;
304
305     //! Callbacks contributing to the conserved energy term
306     std::vector<EnergyContribution> conservedEnergyContributions_;
307     //! Callback to the Parrinello-Rahman box velocities
308     std::function<const rvec*()> parrinelloRahmanBoxVelocities_;
309
310     //! Contains user input mdp options.
311     const t_inputrec* inputrec_;
312     //! Full system topology.
313     const gmx_mtop_t& top_global_;
314     //! Atom parameters for this domain.
315     const MDAtoms* mdAtoms_;
316     //! Energy data structure
317     gmx_enerdata_t* enerd_;
318     //! Kinetic energy data
319     gmx_ekindata_t* ekind_;
320     //! Handles constraints.
321     const Constraints* constr_;
322     //! Handles logging.
323     FILE* fplog_;
324     //! Helper struct for force calculations.
325     t_fcdata* fcd_;
326     //! Notifiers to MD modules
327     const MDModulesNotifiers& mdModulesNotifiers_;
328     //! Global topology groups
329     const SimulationGroups* groups_;
330     //! History of simulation observables.
331     ObservablesHistory* observablesHistory_;
332     //! Whether simulations share the state
333     bool simulationsShareState_;
334 };
335
336 /*! \internal
337  * \ingroup module_modularsimulator
338  * \brief Element for EnergyData
339  *
340  * This member class allows EnergyData to take part in the simulator
341  * loop.
342  *
343  * It subscribes to the trajectory signaller, the energy signaller,
344  * and the logging signaller to know when an energy calculation is
345  * needed and when a non-recording step is enough. The simulator
346  * builder is responsible to place the element in a location at
347  * which a valid energy state is available. The EnergyData::Element is
348  * also a subscriber to the trajectory writer element, as it is
349  * responsible to write energy data to trajectory.
350  */
351 class EnergyData::Element final :
352     public ISimulatorElement,
353     public ITrajectoryWriterClient,
354     public ITrajectorySignallerClient,
355     public IEnergySignallerClient,
356     public ICheckpointHelperClient
357 {
358 public:
359     //! Constructor
360     Element(EnergyData* energyData, bool isMasterRank);
361
362     /*! \brief Register run function for step / time
363      *
364      * This needs to be called when the energies are at a full time step.
365      * Positioning this element is the responsibility of the programmer.
366      *
367      * This is also the place at which the current state becomes the previous
368      * state.
369      *
370      * \param step                 The step number
371      * \param time                 The time
372      * \param registerRunFunction  Function allowing to register a run function
373      */
374     void scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction) override;
375
376     //! No element setup needed
377     void elementSetup() override {}
378
379     //! No element teardown needed
380     void elementTeardown() override {}
381
382     //! ICheckpointHelperClient write checkpoint implementation
383     void saveCheckpointState(std::optional<WriteCheckpointData> checkpointData, const t_commrec* cr) override;
384     //! ICheckpointHelperClient read checkpoint implementation
385     void restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData, const t_commrec* cr) override;
386     //! ICheckpointHelperClient key implementation
387     const std::string& clientID() override;
388
389     /*! \brief Factory method implementation
390      *
391      * \param legacySimulatorData  Pointer allowing access to simulator level data
392      * \param builderHelper  ModularSimulatorAlgorithmBuilder helper object
393      * \param statePropagatorData  Pointer to the \c StatePropagatorData object
394      * \param energyData  Pointer to the \c EnergyData object
395      * \param freeEnergyPerturbationData  Pointer to the \c FreeEnergyPerturbationData object
396      * \param globalCommunicationHelper   Pointer to the \c GlobalCommunicationHelper object
397      * \param observablesReducer          Pointer to the \c ObservablesReducer object
398      *
399      * \return  Pointer to the element to be added. Element needs to have been stored using \c storeElement
400      */
401     static ISimulatorElement* getElementPointerImpl(LegacySimulatorData* legacySimulatorData,
402                                                     ModularSimulatorAlgorithmBuilderHelper* builderHelper,
403                                                     StatePropagatorData*        statePropagatorData,
404                                                     EnergyData*                 energyData,
405                                                     FreeEnergyPerturbationData* freeEnergyPerturbationData,
406                                                     GlobalCommunicationHelper* globalCommunicationHelper,
407                                                     ObservablesReducer*        observablesReducer);
408
409 private:
410     EnergyData* energyData_;
411
412     /*! \brief Setup (needs file pointer)
413      *
414      * ITrajectoryWriterClient implementation.
415      *
416      * Initializes the EnergyOutput object, and does some logging output.
417      *
418      * \param mdoutf  File pointer
419      */
420     void trajectoryWriterSetup(gmx_mdoutf* mdoutf) override;
421     //! No trajectory writer teardown needed
422     void trajectoryWriterTeardown(gmx_mdoutf gmx_unused* outf) override {}
423
424     //! ITrajectoryWriterClient implementation.
425     std::optional<SignallerCallback> registerTrajectorySignallerCallback(TrajectoryEvent event) override;
426     //! ITrajectorySignallerClient implementation
427     std::optional<ITrajectoryWriterCallback> registerTrajectoryWriterCallback(TrajectoryEvent event) override;
428     //! IEnergySignallerClient implementation
429     std::optional<SignallerCallback> registerEnergyCallback(EnergySignallerEvent event) override;
430
431
432     //! CheckpointHelper identifier
433     const std::string identifier_ = "EnergyElement";
434     //! Helper function to read from / write to CheckpointData
435     template<CheckpointDataOperation operation>
436     void doCheckpointData(CheckpointData<operation>* checkpointData);
437
438     //! Whether this is the master rank
439     const bool isMasterRank_;
440     //! The next communicated energy writing step
441     Step energyWritingStep_;
442     //! The next communicated energy calculation step
443     Step energyCalculationStep_;
444     //! The next communicated free energy calculation step
445     Step freeEnergyCalculationStep_;
446 };
447
448 } // namespace gmx
449
450 #endif // GMX_ENERGYELEMENT_MICROSTATE_H