bc0c55dc3081cf32c004721217ee7d9ee1a2b507
[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/math/vectypes.h"
48 #include "gromacs/mdtypes/state.h"
49
50 #include "modularsimulatorinterfaces.h"
51
52 struct gmx_ekindata_t;
53 struct gmx_enerdata_t;
54 struct gmx_mtop_t;
55 struct ObservablesHistory;
56 struct t_fcdata;
57 struct t_inputrec;
58 struct SimulationGroups;
59
60 namespace gmx
61 {
62 enum class StartingBehavior;
63 class Constraints;
64 class EnergyOutput;
65 class FreeEnergyPerturbationData;
66 class GlobalCommunicationHelper;
67 class LegacySimulatorData;
68 class MDAtoms;
69 class ModularSimulatorAlgorithmBuilderHelper;
70 class ParrinelloRahmanBarostat;
71 class StatePropagatorData;
72 class VelocityScalingTemperatureCoupling;
73 struct MdModulesNotifier;
74
75 /*! \internal
76  * \ingroup module_modularsimulator
77  * \brief Data class managing energies
78  *
79  * The EnergyData owns the EnergyObject,
80  * the tensors for the different virials and the pressure as well as
81  * the total dipole vector. It has a member class which is part of the
82  * simulator loop and and is responsible
83  * for saving energy data and writing it to trajectory.
84  *
85  * The EnergyData offers an interface to add virial contributions,
86  * but also allows access to the raw pointers to tensor data, the
87  * dipole vector, and the legacy energy data structures.
88  *
89  * The EnergyData owns an object of type EnergyData::Element,
90  * which takes part in the simulation loop, allowing to record
91  * and output energies during the simulation.
92  */
93 class EnergyData final
94 {
95 public:
96     //! Constructor
97     EnergyData(StatePropagatorData*        statePropagatorData,
98                FreeEnergyPerturbationData* freeEnergyPerturbationData,
99                const gmx_mtop_t&           globalTopology,
100                const t_inputrec*           inputrec,
101                const MDAtoms*              mdAtoms,
102                gmx_enerdata_t*             enerd,
103                gmx_ekindata_t*             ekind,
104                const Constraints*          constr,
105                FILE*                       fplog,
106                t_fcdata*                   fcd,
107                const MdModulesNotifier&    mdModulesNotifier,
108                bool                        isMasterRank,
109                ObservablesHistory*         observablesHistory,
110                StartingBehavior            startingBehavior,
111                bool                        simulationsShareState);
112
113     /*! \brief Final output
114      *
115      * Prints the averages to log. This is called from ModularSimulatorAlgorithm.
116      *
117      * \see ModularSimulatorAlgorithm::teardown
118      */
119     void teardown();
120
121     /*! \brief Add contribution to force virial
122      *
123      * This automatically resets the tensor if the step is higher
124      * than the current step, starting the tensor calculation for
125      * a new step at zero. Otherwise, it adds the new contribution
126      * to the existing virial.
127      */
128     void addToForceVirial(const tensor virial, Step step);
129
130     /*! \brief Add contribution to constraint virial
131      *
132      * This automatically resets the tensor if the step is higher
133      * than the current step, starting the tensor calculation for
134      * a new step at zero. Otherwise, it adds the new contribution
135      * to the existing virial.
136      */
137     void addToConstraintVirial(const tensor virial, Step step);
138
139     /*! \brief Get pointer to force virial tensor
140      *
141      * Allows access to the raw pointer to the tensor.
142      */
143     rvec* forceVirial(Step step);
144
145     /*! \brief Get pointer to constraint virial tensor
146      *
147      * Allows access to the raw pointer to the tensor.
148      */
149     rvec* constraintVirial(Step step);
150
151     /*! \brief Get pointer to total virial tensor
152      *
153      * Allows access to the raw pointer to the tensor.
154      */
155     rvec* totalVirial(Step step);
156
157     /*! \brief Get pointer to pressure tensor
158      *
159      * Allows access to the raw pointer to the tensor.
160      */
161     rvec* pressure(Step step);
162
163     /*! \brief Get pointer to mu_tot
164      *
165      * Allows access to the raw pointer to the dipole vector.
166      */
167     real* muTot();
168
169     /*! \brief Get pointer to energy structure
170      *
171      */
172     gmx_enerdata_t* enerdata();
173
174     /*! \brief Get pointer to kinetic energy structure
175      *
176      */
177     gmx_ekindata_t* ekindata();
178
179     /*! \brief Get pointer to needToSumEkinhOld
180      *
181      */
182     bool* needToSumEkinhOld();
183
184     /*! \brief Whether kinetic energy was read from checkpoint
185      *
186      * This is needed by the compute globals element
187      * TODO: Remove this when moving global reduction to client system (#3421)
188      */
189     [[nodiscard]] bool hasReadEkinFromCheckpoint() const;
190
191     /*! \brief Set velocity scaling temperature coupling
192      *
193      * This allows to set a pointer to a velocity scaling temperature coupling
194      * element used to obtain contributions to the conserved energy.
195      * TODO: This should be made obsolete my a more modular energy element
196      */
197     void setVelocityScalingTemperatureCoupling(const VelocityScalingTemperatureCoupling* velocityScalingTemperatureCoupling);
198
199     /*! \brief set Parrinello-Rahman barostat
200      *
201      * This allows to set a pointer to the Parrinello-Rahman barostat used to
202      * print the box velocities.
203      * TODO: This should be made obsolete my a more modular energy element
204      */
205     void setParrinelloRahamnBarostat(const ParrinelloRahmanBarostat* parrinelloRahmanBarostat);
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(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     //! Pointer to the vrescale thermostat
294     const VelocityScalingTemperatureCoupling* velocityScalingTemperatureCoupling_;
295     //! Pointer to the Parrinello-Rahman barostat
296     const ParrinelloRahmanBarostat* parrinelloRahmanBarostat_;
297     //! Contains user input mdp options.
298     const t_inputrec* inputrec_;
299     //! Full system topology.
300     const gmx_mtop_t& top_global_;
301     //! Atom parameters for this domain.
302     const MDAtoms* mdAtoms_;
303     //! Energy data structure
304     gmx_enerdata_t* enerd_;
305     //! Kinetic energy data
306     gmx_ekindata_t* ekind_;
307     //! Handles constraints.
308     const Constraints* constr_;
309     //! Handles logging.
310     FILE* fplog_;
311     //! Helper struct for force calculations.
312     t_fcdata* fcd_;
313     //! Notification to MD modules
314     const MdModulesNotifier& mdModulesNotifier_;
315     //! Global topology groups
316     const SimulationGroups* groups_;
317     //! History of simulation observables.
318     ObservablesHistory* observablesHistory_;
319     //! Whether simulations share the state
320     bool simulationsShareState_;
321 };
322
323 /*! \internal
324  * \ingroup module_modularsimulator
325  * \brief Element for EnergyData
326  *
327  * This member class allows EnergyData to take part in the simulator
328  * loop.
329  *
330  * It subscribes to the trajectory signaller, the energy signaller,
331  * and the logging signaller to know when an energy calculation is
332  * needed and when a non-recording step is enough. The simulator
333  * builder is responsible to place the element in a location at
334  * which a valid energy state is available. The EnergyData::Element is
335  * also a subscriber to the trajectory writer element, as it is
336  * responsible to write energy data to trajectory.
337  */
338 class EnergyData::Element final :
339     public ISimulatorElement,
340     public ITrajectoryWriterClient,
341     public ITrajectorySignallerClient,
342     public IEnergySignallerClient,
343     public ICheckpointHelperClient
344 {
345 public:
346     //! Constructor
347     Element(EnergyData* energyData, bool isMasterRank);
348
349     /*! \brief Register run function for step / time
350      *
351      * This needs to be called when the energies are at a full time step.
352      * Positioning this element is the responsibility of the programmer.
353      *
354      * This is also the place at which the current state becomes the previous
355      * state.
356      *
357      * \param step                 The step number
358      * \param time                 The time
359      * \param registerRunFunction  Function allowing to register a run function
360      */
361     void scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction) override;
362
363     //! No element setup needed
364     void elementSetup() override {}
365
366     //! No element teardown needed
367     void elementTeardown() override {}
368
369     //! ICheckpointHelperClient write checkpoint implementation
370     void saveCheckpointState(std::optional<WriteCheckpointData> checkpointData, const t_commrec* cr) override;
371     //! ICheckpointHelperClient read checkpoint implementation
372     void restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData, const t_commrec* cr) override;
373     //! ICheckpointHelperClient key implementation
374     const std::string& clientID() override;
375
376     /*! \brief Factory method implementation
377      *
378      * \param legacySimulatorData  Pointer allowing access to simulator level data
379      * \param builderHelper  ModularSimulatorAlgorithmBuilder helper object
380      * \param statePropagatorData  Pointer to the \c StatePropagatorData object
381      * \param energyData  Pointer to the \c EnergyData object
382      * \param freeEnergyPerturbationData  Pointer to the \c FreeEnergyPerturbationData object
383      * \param globalCommunicationHelper  Pointer to the \c GlobalCommunicationHelper object
384      *
385      * \return  Pointer to the element to be added. Element needs to have been stored using \c storeElement
386      */
387     static ISimulatorElement* getElementPointerImpl(LegacySimulatorData* legacySimulatorData,
388                                                     ModularSimulatorAlgorithmBuilderHelper* builderHelper,
389                                                     StatePropagatorData*        statePropagatorData,
390                                                     EnergyData*                 energyData,
391                                                     FreeEnergyPerturbationData* freeEnergyPerturbationData,
392                                                     GlobalCommunicationHelper* globalCommunicationHelper);
393
394 private:
395     EnergyData* energyData_;
396
397     /*! \brief Setup (needs file pointer)
398      *
399      * ITrajectoryWriterClient implementation.
400      *
401      * Initializes the EnergyOutput object, and does some logging output.
402      *
403      * \param mdoutf  File pointer
404      */
405     void trajectoryWriterSetup(gmx_mdoutf* mdoutf) override;
406     //! No trajectory writer teardown needed
407     void trajectoryWriterTeardown(gmx_mdoutf gmx_unused* outf) override {}
408
409     //! ITrajectoryWriterClient implementation.
410     std::optional<SignallerCallback> registerTrajectorySignallerCallback(TrajectoryEvent event) override;
411     //! ITrajectorySignallerClient implementation
412     std::optional<ITrajectoryWriterCallback> registerTrajectoryWriterCallback(TrajectoryEvent event) override;
413     //! IEnergySignallerClient implementation
414     std::optional<SignallerCallback> registerEnergyCallback(EnergySignallerEvent event) override;
415
416
417     //! CheckpointHelper identifier
418     const std::string identifier_ = "EnergyElement";
419     //! Helper function to read from / write to CheckpointData
420     template<CheckpointDataOperation operation>
421     void doCheckpointData(CheckpointData<operation>* checkpointData);
422
423     //! Whether this is the master rank
424     const bool isMasterRank_;
425     //! The next communicated energy writing step
426     Step energyWritingStep_;
427     //! The next communicated energy calculation step
428     Step energyCalculationStep_;
429     //! The next communicated free energy calculation step
430     Step freeEnergyCalculationStep_;
431 };
432
433 } // namespace gmx
434
435 #endif // GMX_ENERGYELEMENT_MICROSTATE_H