Apply clang-format to source tree
[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
51 struct gmx_mdoutf;
52 struct t_commrec;
53 struct t_inputrec;
54 class t_state;
55 struct t_mdatoms;
56
57 namespace gmx
58 {
59 enum class ConstraintVariable;
60 class FreeEnergyPerturbationElement;
61
62 /*! \libinternal
63  * \ingroup module_modularsimulator
64  * \brief StatePropagatorData and associated data
65  *
66  * The `StatePropagatorData` contains a little more than the pure
67  * statistical-physical micro state, namely the positions,
68  * velocities, forces, and box matrix, as well as a backup of
69  * the positions and box of the last time step. While it takes
70  * part in the simulator loop to be able to backup positions /
71  * boxes and save the current state if needed, it's main purpose
72  * is to offer access to its data via getter methods. All elements
73  * reading or writing to this data need a pointer to the
74  * `StatePropagatorData` and need to request their data explicitly. This
75  * will later simplify the understanding of data dependencies
76  * between elements.
77  *
78  * The `StatePropagatorData` takes part in the simulator run, as it might
79  * have to save a valid state at the right moment during the
80  * integration. Placing the StatePropagatorData correctly is for now the
81  * duty of the simulator builder - this might be automatized later
82  * if we have enough meta-data of the variables (i.e., if
83  * `StatePropagatorData` knows at which time the variables currently are,
84  * and can decide when a valid state (full-time step of all
85  * variables) is reached. The `StatePropagatorData` is also a client of
86  * both the trajectory signaller and writer - it will save a
87  * state for later writeout during the simulator step if it
88  * knows that trajectory writing will occur later in the step,
89  * and it knows how to write to file given a file pointer by
90  * the `TrajectoryElement`.
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     public ISimulatorElement,
101     public ITrajectoryWriterClient,
102     public ITrajectorySignallerClient,
103     public ICheckpointHelperClient
104 {
105 public:
106     //! Constructor
107     StatePropagatorData(int                            numAtoms,
108                         FILE*                          fplog,
109                         const t_commrec*               cr,
110                         t_state*                       globalState,
111                         int                            nstxout,
112                         int                            nstvout,
113                         int                            nstfout,
114                         int                            nstxout_compressed,
115                         bool                           useGPU,
116                         FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
117                         const t_inputrec*              inputrec,
118                         const t_mdatoms*               mdatoms);
119
120     // Allow access to state
121     //! Get write access to position vector
122     ArrayRefWithPadding<RVec> positionsView();
123     //! Get read access to position vector
124     ArrayRefWithPadding<const RVec> constPositionsView() const;
125     //! Get write access to previous position vector
126     ArrayRefWithPadding<RVec> previousPositionsView();
127     //! Get read access to previous position vector
128     ArrayRefWithPadding<const RVec> constPreviousPositionsView() const;
129     //! Get write access to velocity vector
130     ArrayRefWithPadding<RVec> velocitiesView();
131     //! Get read access to velocity vector
132     ArrayRefWithPadding<const RVec> constVelocitiesView() const;
133     //! Get write access to force vector
134     ArrayRefWithPadding<RVec> forcesView();
135     //! Get read access to force vector
136     ArrayRefWithPadding<const RVec> constForcesView() const;
137     //! Get pointer to box
138     rvec* box();
139     //! Get const pointer to box
140     const rvec* constBox();
141     //! Get pointer to previous box
142     rvec* previousBox();
143     //! Get const pointer to previous box
144     const rvec* constPreviousBox();
145     //! Get the local number of atoms
146     int localNumAtoms();
147
148     /*! \brief Register run function for step / time
149      *
150      * This needs to be called during the integration part of the simulator,
151      * at the moment at which the state is at a full time step. Positioning
152      * this element is the responsibility of the programmer writing the
153      * integration algorithm! If the current step is a trajectory writing
154      * step, StatePropagatorData will save a backup for later writeout.
155      *
156      * This is also the place at which the current state becomes the previous
157      * state.
158      *
159      * @param step                 The step number
160      * @param time                 The time
161      * @param registerRunFunction  Function allowing to register a run function
162      */
163     void scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction) override;
164
165     /*! \brief Backup starting velocities
166      *
167      * This is only needed for vv, where the first (velocity) half step is only
168      * used to compute the constraint virial, but the velocities need to be reset
169      * after.
170      * TODO: There must be a more elegant solution to this!
171      */
172     void elementSetup() override;
173
174     //! No element teardown needed
175     void elementTeardown() override {}
176
177     //! @cond
178     // (doxygen doesn't like these)
179     // Classes which need access to legacy state
180     friend class DomDecHelper;
181     //! @endcond
182
183 private:
184     //! The total number of atoms in the system
185     int totalNumAtoms_;
186     //! The position writeout frequency
187     int nstxout_;
188     //! The velocity writeout frequency
189     int nstvout_;
190     //! The force writeout frequency
191     int nstfout_;
192     //! The compressed position writeout frequency
193     int nstxout_compressed_;
194
195     //! The local number of atoms
196     int localNAtoms_;
197     //! The position vector
198     PaddedHostVector<RVec> x_;
199     //! The position vector of the previous step
200     PaddedHostVector<RVec> previousX_;
201     //! The velocity vector
202     PaddedHostVector<RVec> v_;
203     //! The force vector
204     PaddedHostVector<RVec> f_;
205     //! The box matrix
206     matrix box_;
207     //! The box matrix of the previous step
208     matrix previousBox_;
209     //! The DD partitioning count for legacy t_state compatibility
210     int ddpCount_;
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     // Access to legacy state
218     //! Get a deep copy of the current state in legacy format
219     std::unique_ptr<t_state> localState();
220     //! Update the current state with a state in legacy format
221     void setLocalState(std::unique_ptr<t_state> state);
222     //! Get a pointer to the global state
223     t_state* globalState();
224     //! Get a force pointer
225     PaddedHostVector<gmx::RVec>* forcePointer();
226
227     //! Pointer to keep a backup of the state for later writeout
228     std::unique_ptr<t_state> localStateBackup_;
229     //! Step at which next writeout occurs
230     Step writeOutStep_;
231     //! Backup current state
232     void saveState();
233
234     //! ITrajectorySignallerClient implementation
235     SignallerCallbackPtr registerTrajectorySignallerCallback(TrajectoryEvent event) override;
236
237     //! ITrajectoryWriterClient implementation
238     ITrajectoryWriterCallbackPtr registerTrajectoryWriterCallback(TrajectoryEvent event) override;
239
240     //! ICheckpointHelperClient implementation
241     void writeCheckpoint(t_state* localState, t_state* globalState) override;
242
243     //! Callback writing the state to file
244     void write(gmx_mdoutf* outf, Step step, Time time);
245
246     //! Whether we're doing VV and need to reset velocities after the first half step
247     bool vvResetVelocities_;
248     //! Velocities backup for VV
249     PaddedHostVector<RVec> velocityBackup_;
250     //! Function resetting the velocities
251     void resetVelocities();
252
253     //! Pointer to the free energy perturbation element (for trajectory writing only)
254     FreeEnergyPerturbationElement* freeEnergyPerturbationElement_;
255
256     // Access to ISimulator data
257     //! Handles logging.
258     FILE* fplog_;
259     //! Handles communication.
260     const t_commrec* cr_;
261     //! Full simulation state (only non-nullptr on master rank).
262     t_state* globalState_;
263
264     //! No trajectory writer setup needed
265     void trajectoryWriterSetup(gmx_mdoutf gmx_unused* outf) override {}
266     //! No trajectory writer teardown needed
267     void trajectoryWriterTeardown(gmx_mdoutf gmx_unused* outf) override {}
268 };
269
270 } // namespace gmx
271
272 #endif // GMX_MODULARSIMULATOR_STATEPROPAGATORDATA_H