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