Domain decomposition and PME load balancing for modular simulator
[gromacs.git] / src / gromacs / modularsimulator / statepropagatordata.h
blob20dbe13cb5a607c012252b2300219fa92331ecec
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
35 /*! \libinternal
36 * \brief Declares the state for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
42 #ifndef GMX_MODULARSIMULATOR_STATEPROPAGATORDATA_H
43 #define GMX_MODULARSIMULATOR_STATEPROPAGATORDATA_H
45 #include "gromacs/gpu_utils/hostallocator.h"
46 #include "gromacs/math/paddedvector.h"
47 #include "gromacs/math/vectypes.h"
49 #include "modularsimulatorinterfaces.h"
51 struct gmx_mdoutf;
52 struct t_commrec;
53 struct t_inputrec;
54 class t_state;
55 struct t_mdatoms;
57 namespace gmx
59 enum class ConstraintVariable;
61 //! \addtogroup module_modularsimulator
62 //! \{
64 /*! \libinternal
65 * \brief StatePropagatorData and associated data
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.
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`.
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.
100 class StatePropagatorData final :
101 public ISimulatorElement,
102 public ITrajectoryWriterClient,
103 public ITrajectorySignallerClient
105 public:
106 //! Constructor
107 StatePropagatorData(
108 int numAtoms,
109 FILE *fplog,
110 const t_commrec *cr,
111 t_state *globalState,
112 int nstxout,
113 int nstvout,
114 int nstfout,
115 int nstxout_compressed,
116 bool useGPU,
117 const t_inputrec *inputrec,
118 const t_mdatoms *mdatoms);
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 pointer to previous box
140 rvec *previousBox();
141 //! Get the local number of atoms
142 int localNumAtoms();
144 /*! \brief Register run function for step / time
146 * This needs to be called during the integration part of the simulator,
147 * at the moment at which the state is at a full time step. Positioning
148 * this element is the responsibility of the programmer writing the
149 * integration algorithm! If the current step is a trajectory writing
150 * step, StatePropagatorData will save a backup for later writeout.
152 * This is also the place at which the current state becomes the previous
153 * state.
155 * @param step The step number
156 * @param time The time
157 * @param registerRunFunction Function allowing to register a run function
159 void scheduleTask(
160 Step step, Time time,
161 const RegisterRunFunctionPtr &registerRunFunction) override;
163 /*! \brief Backup starting velocities
165 * This is only needed for vv, where the first (velocity) half step is only
166 * used to compute the constraint virial, but the velocities need to be reset
167 * after.
168 * TODO: There must be a more elegant solution to this!
170 void elementSetup() override;
172 //! No element teardown needed
173 void elementTeardown() override {}
175 //! @cond
176 // (doxygen doesn't like these)
177 // Classes which need access to legacy state
178 friend class DomDecHelper;
179 friend class PmeLoadBalanceHelper;
180 //! @endcond
182 private:
183 //! The total number of atoms in the system
184 int totalNumAtoms_;
185 //! The position writeout frequency
186 int nstxout_;
187 //! The velocity writeout frequency
188 int nstvout_;
189 //! The force writeout frequency
190 int nstfout_;
191 //! The compressed position writeout frequency
192 int nstxout_compressed_;
194 //! The local number of atoms
195 int localNAtoms_;
196 //! The position vector
197 PaddedHostVector<RVec> x_;
198 //! The position vector of the previous step
199 PaddedVector<RVec> previousX_;
200 //! The velocity vector
201 PaddedVector<RVec> v_;
202 //! The force vector
203 PaddedVector<RVec> f_;
204 //! The box matrix
205 matrix box_;
206 //! The box matrix of the previous step
207 matrix previousBox_;
208 //! Bit-flags for legacy t_state compatibility
209 unsigned int flags_;
210 //! The DD partitioning count for legacy t_state compatibility
211 int ddpCount_;
213 //! Move x_ to previousX_
214 void copyPosition();
215 //! OMP helper to move x_ to previousX_
216 void copyPosition(int start, int end);
218 // Access to legacy state
219 //! Get a deep copy of the current state in legacy format
220 std::unique_ptr<t_state> localState();
221 //! Update the current state with a state in legacy format
222 void setLocalState(std::unique_ptr<t_state> state);
223 //! Get a pointer to the global state
224 t_state *globalState();
225 //! Get a force pointer
226 PaddedVector<gmx::RVec> *forcePointer();
228 //! Pointer to keep a backup of the state for later writeout
229 std::unique_ptr<t_state> localStateBackup_;
230 //! Step at which next writeout occurs
231 Step writeOutStep_;
232 //! Backup current state
233 void saveState();
235 //! ITrajectorySignallerClient implementation
236 SignallerCallbackPtr
237 registerTrajectorySignallerCallback(TrajectoryEvent event) override;
239 //! ITrajectoryWriterClient implementation
240 ITrajectoryWriterCallbackPtr
241 registerTrajectoryWriterCallback(TrajectoryEvent event) override;
243 //! Callback writing the state to file
244 void write(gmx_mdoutf *outf, Step step, Time time);
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 PaddedVector<RVec> velocityBackup_;
250 //! Function resetting the velocities
251 void resetVelocities();
253 // Access to ISimulator data
254 //! Handles logging.
255 FILE *fplog_;
256 //! Handles communication.
257 const t_commrec *cr_;
258 //! Full simulation state (only non-nullptr on master rank).
259 t_state *globalState_;
261 //! No trajectory writer setup needed
262 void trajectoryWriterSetup(gmx_mdoutf gmx_unused *outf) override {}
263 //! No trajectory writer teardown needed
264 void trajectoryWriterTeardown(gmx_mdoutf gmx_unused *outf) override {}
267 //! /}
268 } // namespace gmx
270 #endif // GMX_MODULARSIMULATOR_STATEPROPAGATORDATA_H