78007403d50b0f083d3e268baba1f548329fd84d
[alexxy/gromacs.git] / src / gromacs / modularsimulator / statepropagatordata.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 state 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_MODULARSIMULATOR_STATEPROPAGATORDATA_H
45 #define GMX_MODULARSIMULATOR_STATEPROPAGATORDATA_H
46
47 #include "gromacs/gpu_utils/hostallocator.h"
48 #include "gromacs/math/paddedvector.h"
49 #include "gromacs/math/vectypes.h"
50 #include "gromacs/mdtypes/checkpointdata.h"
51 #include "gromacs/mdtypes/forcebuffers.h"
52 #include "gromacs/utility/keyvaluetree.h"
53
54 #include "modularsimulatorinterfaces.h"
55 #include "topologyholder.h"
56
57 struct gmx_mdoutf;
58 enum class PbcType : int;
59 struct t_commrec;
60 struct t_inputrec;
61 class t_state;
62 struct t_mdatoms;
63 struct t_trxframe;
64
65 namespace gmx
66 {
67 enum class CheckpointDataOperation;
68 enum class ConstraintVariable;
69 class EnergyData;
70 class FreeEnergyPerturbationData;
71 class GlobalCommunicationHelper;
72 class LegacySimulatorData;
73 class ModularSimulatorAlgorithmBuilderHelper;
74
75 /*! \internal
76  * \ingroup module_modularsimulator
77  * \brief StatePropagatorData and associated data
78  *
79  * The `StatePropagatorData` contains a little more than the pure
80  * statistical-physical micro state, namely the positions,
81  * velocities, forces, and box matrix, as well as a backup of
82  * the positions and box of the last time step. While it takes
83  * part in the simulator loop via its member class `Element`
84  * to be able to backup positions /
85  * boxes and save the current state if needed, it's main purpose
86  * is to offer access to its data via getter methods. All elements
87  * reading or writing to this data need a pointer to the
88  * `StatePropagatorData` and need to request their data explicitly. This
89  * will later simplify the understanding of data dependencies
90  * between elements.
91  *
92  * Note that the `StatePropagatorData` can be converted to and from the
93  * legacy `t_state` object. This is useful when dealing with
94  * functionality which has not yet been adapted to use the new
95  * data approach - of the elements currently implemented, only
96  * domain decomposition, PME load balancing, and the initial
97  * constraining are using this.
98  */
99 class StatePropagatorData final
100 {
101 public:
102     //! Constructor
103     StatePropagatorData(int                numAtoms,
104                         FILE*              fplog,
105                         const t_commrec*   cr,
106                         t_state*           globalState,
107                         bool               useGPU,
108                         bool               canMoleculesBeDistributedOverPBC,
109                         bool               writeFinalConfiguration,
110                         const std::string& finalConfigurationFilename,
111                         const t_inputrec*  inputrec,
112                         const t_mdatoms*   mdatoms,
113                         const gmx_mtop_t&  globalTop);
114
115     //! Destructor (allows forward declaration of internal type)
116     ~StatePropagatorData();
117
118     // Allow access to state
119     //! Get write access to position vector
120     ArrayRefWithPadding<RVec> positionsView();
121     //! Get read access to position vector
122     ArrayRefWithPadding<const RVec> constPositionsView() const;
123     //! Get write access to previous position vector
124     ArrayRefWithPadding<RVec> previousPositionsView();
125     //! Get read access to previous position vector
126     ArrayRefWithPadding<const RVec> constPreviousPositionsView() const;
127     //! Get write access to velocity vector
128     ArrayRefWithPadding<RVec> velocitiesView();
129     //! Get read access to velocity vector
130     ArrayRefWithPadding<const RVec> constVelocitiesView() const;
131     //! Get write access to force vector
132     ForceBuffersView& forcesView();
133     //! Get read access to force vector
134     const ForceBuffersView& constForcesView() const;
135     //! Get pointer to box
136     rvec* box();
137     //! Get const pointer to box
138     const rvec* constBox() const;
139     //! Get pointer to previous box
140     rvec* previousBox();
141     //! Get const pointer to previous box
142     const rvec* constPreviousBox() const;
143     //! Get the local number of atoms
144     int localNumAtoms() const;
145     //! Get the total number of atoms
146     int totalNumAtoms() const;
147
148     //! The element taking part in the simulator loop
149     class Element;
150     //! Get pointer to element (whose lifetime is managed by this)
151     Element* element();
152     //! Initial set up for the associated element
153     void setup();
154
155     //! Update the reference temperature
156     void updateReferenceTemperature(ArrayRef<const real>                temperatures,
157                                     ReferenceTemperatureChangeAlgorithm algorithm);
158     //! Helper class handling reference temperature change
159     class ReferenceTemperatureHelper;
160
161     //! Read everything that can be stored in t_trxframe from a checkpoint file
162     static void readCheckpointToTrxFrame(t_trxframe* trxFrame, ReadCheckpointData readCheckpointData);
163     //! CheckpointHelper identifier
164     static const std::string& checkpointID();
165
166     //! \cond
167     // (doxygen doesn't like these)
168     // Classes which need access to legacy state
169     friend class DomDecHelper;
170     //! \endcond
171
172 private:
173     //! Default constructor - only used internally
174     StatePropagatorData() = default;
175     //! The total number of atoms in the system
176     int totalNumAtoms_;
177     //! The local number of atoms
178     int localNAtoms_;
179     //! The position vector
180     PaddedHostVector<RVec> x_;
181     //! The position vector of the previous step
182     PaddedHostVector<RVec> previousX_;
183     //! The velocity vector
184     PaddedHostVector<RVec> v_;
185     //! The force vector
186     ForceBuffers f_;
187     //! The box matrix
188     matrix box_;
189     //! The box matrix of the previous step
190     matrix previousBox_;
191     //! The DD partitioning count
192     int ddpCount_;
193     //! The DD partitioning count for index_gl
194     int ddpCountCgGl_;
195     //! The global cg number of the local cgs
196     std::vector<int> cgGl_;
197
198     //! The global position vector
199     PaddedHostVector<RVec> xGlobal_;
200     //! The global position vector of the previous step
201     PaddedHostVector<RVec> previousXGlobal_;
202     //! The global velocity vector
203     PaddedHostVector<RVec> vGlobal_;
204     //! The global force vector
205     PaddedHostVector<RVec> fGlobal_;
206
207     //! The element
208     std::unique_ptr<Element> element_;
209     //! Instance of reference temperature helper
210     std::unique_ptr<ReferenceTemperatureHelper> referenceTemperatureHelper_;
211
212     //! Move x_ to previousX_
213     void copyPosition();
214     //! OMP helper to move x_ to previousX_
215     void copyPosition(int start, int end);
216
217     //! Helper function to read from / write to CheckpointData
218     template<CheckpointDataOperation operation>
219     void doCheckpointData(CheckpointData<operation>* checkpointData);
220
221     // Access to legacy state
222     //! Get a deep copy of the current state in legacy format
223     std::unique_ptr<t_state> localState();
224     //! Update the current state with a state in legacy format
225     void setLocalState(std::unique_ptr<t_state> state);
226     //! Get a pointer to the global state
227     t_state* globalState();
228     //! Get a force pointer
229     ForceBuffers* forcePointer();
230
231     //! Whether we're doing VV and need to reset velocities after the first half step
232     bool vvResetVelocities_;
233     //! Velocities backup for VV
234     PaddedHostVector<RVec> velocityBackup_;
235     //! Function resetting the velocities
236     void resetVelocities();
237
238     //! Whether planned total number of steps was reached (used for final output only)
239     bool isRegularSimulationEnd_;
240     //! The signalled last step (used for final output only)
241     Step lastStep_;
242
243     // Access to ISimulator data
244     //! Full simulation state (only non-nullptr on master rank).
245     t_state* globalState_;
246 };
247
248 /*! \internal
249  * \ingroup module_modularsimulator
250  * \brief Element for StatePropagatorData
251  *
252  * The `StatePropagatorData::Element` takes part in the simulator run, as it might
253  * have to save a valid state at the right moment during the
254  * integration. Placing the StatePropagatorData::Element correctly is the
255  * duty of the simulator builder - this might be automatized later
256  * if we have enough meta-data of the variables (i.e., if
257  * `StatePropagatorData` knows at which time the variables currently are,
258  * and can decide when a valid state (full-time step of all
259  * variables) is reached. The `StatePropagatorData::Element` is also a client of
260  * both the trajectory signaller and writer - it will save a
261  * state for later writeout during the simulator step if it
262  * knows that trajectory writing will occur later in the step,
263  * and it knows how to write to file given a file pointer by
264  * the `TrajectoryElement`. It is also responsible to store
265  * the state for checkpointing.
266  *
267  */
268 class StatePropagatorData::Element final :
269     public ISimulatorElement,
270     public ITrajectoryWriterClient,
271     public ITrajectorySignallerClient,
272     public ICheckpointHelperClient,
273     public ILastStepSignallerClient
274 {
275 public:
276     //! Constructor
277     Element(StatePropagatorData* statePropagatorData,
278             FILE*                fplog,
279             const t_commrec*     cr,
280             int                  nstxout,
281             int                  nstvout,
282             int                  nstfout,
283             int                  nstxout_compressed,
284             bool                 canMoleculesBeDistributedOverPBC,
285             bool                 writeFinalConfiguration,
286             std::string          finalConfigurationFilename,
287             const t_inputrec*    inputrec,
288             const gmx_mtop_t&    globalTop);
289
290     /*! \brief Register run function for step / time
291      *
292      * This needs to be called during the integration part of the simulator,
293      * at the moment at which the state is at a full time step. Positioning
294      * this element is the responsibility of the programmer writing the
295      * integration algorithm! If the current step is a trajectory writing
296      * step, StatePropagatorData will save a backup for later writeout.
297      *
298      * This is also the place at which the current state becomes the previous
299      * state.
300      *
301      * \param step                 The step number
302      * \param time                 The time
303      * \param registerRunFunction  Function allowing to register a run function
304      */
305     void scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction) override;
306
307     /*! \brief Backup starting velocities
308      *
309      * This is only needed for vv, where the first (velocity) half step is only
310      * used to compute the constraint virial, but the velocities need to be reset
311      * after.
312      * TODO: There must be a more elegant solution to this!
313      */
314     void elementSetup() override;
315
316     //! No element teardown needed
317     void elementTeardown() override {}
318
319     //! Set free energy data
320     void setFreeEnergyPerturbationData(FreeEnergyPerturbationData* freeEnergyPerturbationData);
321
322     //! ICheckpointHelperClient write checkpoint implementation
323     void saveCheckpointState(std::optional<WriteCheckpointData> checkpointData, const t_commrec* cr) override;
324     //! ICheckpointHelperClient read checkpoint implementation
325     void restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData, const t_commrec* cr) override;
326     //! ICheckpointHelperClient key implementation
327     const std::string& clientID() override;
328
329     /*! \brief Factory method implementation
330      *
331      * \param legacySimulatorData  Pointer allowing access to simulator level data
332      * \param builderHelper  ModularSimulatorAlgorithmBuilder helper object
333      * \param statePropagatorData  Pointer to the \c StatePropagatorData object
334      * \param energyData  Pointer to the \c EnergyData object
335      * \param freeEnergyPerturbationData  Pointer to the \c FreeEnergyPerturbationData object
336      * \param globalCommunicationHelper  Pointer to the \c GlobalCommunicationHelper object
337      *
338      * \return  Pointer to the element to be added. Element needs to have been stored using \c storeElement
339      */
340     static ISimulatorElement* getElementPointerImpl(LegacySimulatorData* legacySimulatorData,
341                                                     ModularSimulatorAlgorithmBuilderHelper* builderHelper,
342                                                     StatePropagatorData*        statePropagatorData,
343                                                     EnergyData*                 energyData,
344                                                     FreeEnergyPerturbationData* freeEnergyPerturbationData,
345                                                     GlobalCommunicationHelper* globalCommunicationHelper);
346
347 private:
348     //! Pointer to the associated StatePropagatorData
349     StatePropagatorData* statePropagatorData_;
350
351     //! The position writeout frequency
352     const int nstxout_;
353     //! The velocity writeout frequency
354     const int nstvout_;
355     //! The force writeout frequency
356     const int nstfout_;
357     //! The compressed position writeout frequency
358     const int nstxout_compressed_;
359
360     //! Pointer to keep a backup of the state for later writeout
361     std::unique_ptr<t_state> localStateBackup_;
362     //! Step at which next writeout occurs
363     Step writeOutStep_;
364     //! Backup current state
365     void saveState();
366
367     //! ITrajectorySignallerClient implementation
368     std::optional<SignallerCallback> registerTrajectorySignallerCallback(TrajectoryEvent event) override;
369
370     //! ITrajectoryWriterClient implementation
371     std::optional<ITrajectoryWriterCallback> registerTrajectoryWriterCallback(TrajectoryEvent event) override;
372
373     //! ILastStepSignallerClient implementation (used for final output only)
374     std::optional<SignallerCallback> registerLastStepCallback() override;
375
376     //! Callback writing the state to file
377     void write(gmx_mdoutf* outf, Step step, Time time);
378
379     // TODO: Clarify relationship to data objects and find a more robust alternative to raw pointers (#3583)
380     //! Pointer to the free energy perturbation data (for trajectory writing only)
381     FreeEnergyPerturbationData* freeEnergyPerturbationData_;
382
383     //! No trajectory writer setup needed
384     void trajectoryWriterSetup(gmx_mdoutf gmx_unused* outf) override {}
385     //! Trajectory writer teardown - write final coordinates
386     void trajectoryWriterTeardown(gmx_mdoutf* outf) override;
387     //! A dummy CheckpointData - remove when we stop using the legacy trajectory writing function
388     WriteCheckpointDataHolder dummyCheckpointDataHolder_;
389
390     //! Whether planned total number of steps was reached (used for final output only)
391     bool isRegularSimulationEnd_;
392     //! The signalled last step (used for final output only)
393     Step lastStep_;
394     //! Whether system can have molecules distributed over PBC boundaries (used for final output only)
395     const bool canMoleculesBeDistributedOverPBC_;
396     //! Whether system has molecules self-interacting through PBC (used for final output only)
397     const bool systemHasPeriodicMolecules_;
398     //! The PBC type (used for final output only)
399     const PbcType pbcType_;
400     //! The (planned) last step - determines whether final configuration is written (used for final output only)
401     const Step lastPlannedStep_;
402     //! Whether final configuration was chosen in mdrun options (used for final output only)
403     const bool writeFinalConfiguration_;
404     //! The filename of the final configuration file (used for final output only)
405     const std::string finalConfigurationFilename_;
406
407     // Access to ISimulator data
408     //! Handles logging.
409     FILE* fplog_;
410     //! Handles communication.
411     const t_commrec* cr_;
412     //! Full system topology.
413     const gmx_mtop_t& top_global_;
414 };
415
416 } // namespace gmx
417
418 #endif // GMX_MODULARSIMULATOR_STATEPROPAGATORDATA_H