Change the behavior of the GPU update UI
[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, 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 /*! \libinternal \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
42 #ifndef GMX_MODULARSIMULATOR_STATEPROPAGATORDATA_H
43 #define GMX_MODULARSIMULATOR_STATEPROPAGATORDATA_H
44
45 #include "gromacs/gpu_utils/hostallocator.h"
46 #include "gromacs/math/paddedvector.h"
47 #include "gromacs/math/vectypes.h"
48
49 #include "modularsimulatorinterfaces.h"
50 #include "topologyholder.h"
51
52 struct gmx_mdoutf;
53 struct t_commrec;
54 struct t_inputrec;
55 class t_state;
56 struct t_mdatoms;
57
58 namespace gmx
59 {
60 enum class ConstraintVariable;
61 class FreeEnergyPerturbationElement;
62
63 /*! \libinternal
64  * \ingroup module_modularsimulator
65  * \brief StatePropagatorData and associated data
66  *
67  * The `StatePropagatorData` contains a little more than the pure
68  * statistical-physical micro state, namely the positions,
69  * velocities, forces, and box matrix, as well as a backup of
70  * the positions and box of the last time step. While it takes
71  * part in the simulator loop to be able to backup positions /
72  * boxes and save the current state if needed, it's main purpose
73  * is to offer access to its data via getter methods. All elements
74  * reading or writing to this data need a pointer to the
75  * `StatePropagatorData` and need to request their data explicitly. This
76  * will later simplify the understanding of data dependencies
77  * between elements.
78  *
79  * The `StatePropagatorData` takes part in the simulator run, as it might
80  * have to save a valid state at the right moment during the
81  * integration. Placing the StatePropagatorData correctly is for now the
82  * duty of the simulator builder - this might be automatized later
83  * if we have enough meta-data of the variables (i.e., if
84  * `StatePropagatorData` knows at which time the variables currently are,
85  * and can decide when a valid state (full-time step of all
86  * variables) is reached. The `StatePropagatorData` is also a client of
87  * both the trajectory signaller and writer - it will save a
88  * state for later writeout during the simulator step if it
89  * knows that trajectory writing will occur later in the step,
90  * and it knows how to write to file given a file pointer by
91  * the `TrajectoryElement`.
92  *
93  * Note that the `StatePropagatorData` can be converted to and from the
94  * legacy `t_state` object. This is useful when dealing with
95  * functionality which has not yet been adapted to use the new
96  * data approach - of the elements currently implemented, only
97  * domain decomposition, PME load balancing, and the initial
98  * constraining are using this.
99  */
100 class StatePropagatorData final :
101     public ISimulatorElement,
102     public ITrajectoryWriterClient,
103     public ITrajectorySignallerClient,
104     public ICheckpointHelperClient,
105     public ILastStepSignallerClient
106 {
107 public:
108     //! Constructor
109     StatePropagatorData(int                            numAtoms,
110                         FILE*                          fplog,
111                         const t_commrec*               cr,
112                         t_state*                       globalState,
113                         int                            nstxout,
114                         int                            nstvout,
115                         int                            nstfout,
116                         int                            nstxout_compressed,
117                         bool                           useGPU,
118                         FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
119                         const TopologyHolder*          topologyHolder,
120                         bool                           canMoleculesBeDistributedOverPBC,
121                         bool                           writeFinalConfiguration,
122                         std::string                    finalConfigurationFilename,
123                         const t_inputrec*              inputrec,
124                         const t_mdatoms*               mdatoms);
125
126     // Allow access to state
127     //! Get write access to position vector
128     ArrayRefWithPadding<RVec> positionsView();
129     //! Get read access to position vector
130     ArrayRefWithPadding<const RVec> constPositionsView() const;
131     //! Get write access to previous position vector
132     ArrayRefWithPadding<RVec> previousPositionsView();
133     //! Get read access to previous position vector
134     ArrayRefWithPadding<const RVec> constPreviousPositionsView() const;
135     //! Get write access to velocity vector
136     ArrayRefWithPadding<RVec> velocitiesView();
137     //! Get read access to velocity vector
138     ArrayRefWithPadding<const RVec> constVelocitiesView() const;
139     //! Get write access to force vector
140     ArrayRefWithPadding<RVec> forcesView();
141     //! Get read access to force vector
142     ArrayRefWithPadding<const RVec> constForcesView() const;
143     //! Get pointer to box
144     rvec* box();
145     //! Get const pointer to box
146     const rvec* constBox();
147     //! Get pointer to previous box
148     rvec* previousBox();
149     //! Get const pointer to previous box
150     const rvec* constPreviousBox();
151     //! Get the local number of atoms
152     int localNumAtoms();
153
154     /*! \brief Register run function for step / time
155      *
156      * This needs to be called during the integration part of the simulator,
157      * at the moment at which the state is at a full time step. Positioning
158      * this element is the responsibility of the programmer writing the
159      * integration algorithm! If the current step is a trajectory writing
160      * step, StatePropagatorData will save a backup for later writeout.
161      *
162      * This is also the place at which the current state becomes the previous
163      * state.
164      *
165      * @param step                 The step number
166      * @param time                 The time
167      * @param registerRunFunction  Function allowing to register a run function
168      */
169     void scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction) override;
170
171     /*! \brief Backup starting velocities
172      *
173      * This is only needed for vv, where the first (velocity) half step is only
174      * used to compute the constraint virial, but the velocities need to be reset
175      * after.
176      * TODO: There must be a more elegant solution to this!
177      */
178     void elementSetup() override;
179
180     //! No element teardown needed
181     void elementTeardown() override {}
182
183     //! @cond
184     // (doxygen doesn't like these)
185     // Classes which need access to legacy state
186     friend class DomDecHelper;
187     //! @endcond
188
189 private:
190     //! The total number of atoms in the system
191     int totalNumAtoms_;
192     //! The position writeout frequency
193     int nstxout_;
194     //! The velocity writeout frequency
195     int nstvout_;
196     //! The force writeout frequency
197     int nstfout_;
198     //! The compressed position writeout frequency
199     int nstxout_compressed_;
200
201     //! The local number of atoms
202     int localNAtoms_;
203     //! The position vector
204     PaddedHostVector<RVec> x_;
205     //! The position vector of the previous step
206     PaddedHostVector<RVec> previousX_;
207     //! The velocity vector
208     PaddedHostVector<RVec> v_;
209     //! The force vector
210     PaddedHostVector<RVec> f_;
211     //! The box matrix
212     matrix box_;
213     //! The box matrix of the previous step
214     matrix previousBox_;
215     //! The DD partitioning count for legacy t_state compatibility
216     int ddpCount_;
217
218     //! Move x_ to previousX_
219     void copyPosition();
220     //! OMP helper to move x_ to previousX_
221     void copyPosition(int start, int end);
222
223     // Access to legacy state
224     //! Get a deep copy of the current state in legacy format
225     std::unique_ptr<t_state> localState();
226     //! Update the current state with a state in legacy format
227     void setLocalState(std::unique_ptr<t_state> state);
228     //! Get a pointer to the global state
229     t_state* globalState();
230     //! Get a force pointer
231     PaddedHostVector<gmx::RVec>* forcePointer();
232
233     //! Pointer to keep a backup of the state for later writeout
234     std::unique_ptr<t_state> localStateBackup_;
235     //! Step at which next writeout occurs
236     Step writeOutStep_;
237     //! Backup current state
238     void saveState();
239
240     //! ITrajectorySignallerClient implementation
241     SignallerCallbackPtr registerTrajectorySignallerCallback(TrajectoryEvent event) override;
242
243     //! ITrajectoryWriterClient implementation
244     ITrajectoryWriterCallbackPtr registerTrajectoryWriterCallback(TrajectoryEvent event) override;
245
246     //! ICheckpointHelperClient implementation
247     void writeCheckpoint(t_state* localState, t_state* globalState) override;
248
249     //! ILastStepSignallerClient implementation (used for final output only)
250     SignallerCallbackPtr registerLastStepCallback() override;
251
252     //! Callback writing the state to file
253     void write(gmx_mdoutf* outf, Step step, Time time);
254
255     //! Whether we're doing VV and need to reset velocities after the first half step
256     bool vvResetVelocities_;
257     //! Velocities backup for VV
258     PaddedHostVector<RVec> velocityBackup_;
259     //! Function resetting the velocities
260     void resetVelocities();
261
262     //! Pointer to the free energy perturbation element (for trajectory writing only)
263     FreeEnergyPerturbationElement* freeEnergyPerturbationElement_;
264
265     //! Whether planned total number of steps was reached (used for final output only)
266     bool isRegularSimulationEnd_;
267     //! The signalled last step (used for final output only)
268     Step lastStep_;
269
270     //! Whether system can have molecules distributed over PBC boundaries (used for final output only)
271     const bool canMoleculesBeDistributedOverPBC_;
272     //! Whether system has molecules self-interacting through PBC (used for final output only)
273     const bool systemHasPeriodicMolecules_;
274     //! The PBC type (used for final output only)
275     const int pbcType_;
276     //! Pointer to the topology (used for final output only)
277     const TopologyHolder* topologyHolder_;
278     //! The (planned) last step - determines whether final configuration is written (used for final output only)
279     const Step lastPlannedStep_;
280     //! Whether final configuration was chosen in mdrun options (used for final output only)
281     const bool writeFinalConfiguration_;
282     //! The filename of the final configuration file (used for final output only)
283     const std::string finalConfigurationFilename_;
284
285     // Access to ISimulator data
286     //! Handles logging.
287     FILE* fplog_;
288     //! Handles communication.
289     const t_commrec* cr_;
290     //! Full simulation state (only non-nullptr on master rank).
291     t_state* globalState_;
292
293     //! No trajectory writer setup needed
294     void trajectoryWriterSetup(gmx_mdoutf gmx_unused* outf) override {}
295     //! Trajectory writer teardown - write final coordinates
296     void trajectoryWriterTeardown(gmx_mdoutf* outf) override;
297 };
298
299 } // namespace gmx
300
301 #endif // GMX_MODULARSIMULATOR_STATEPROPAGATORDATA_H