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