Remove unused thole polarization rfac parameter
[alexxy/gromacs.git] / src / gromacs / fileio / checkpoint.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38
39 #ifndef GMX_FILEIO_CHECKPOINT_H
40 #define GMX_FILEIO_CHECKPOINT_H
41
42 #include <cstdio>
43
44 #include <vector>
45
46 #include "gromacs/compat/pointers.h"
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/utility/basedefinitions.h"
49 #include "gromacs/utility/gmxmpi.h"
50 #include "gromacs/utility/keyvaluetreebuilder.h"
51
52 class energyhistory_t;
53 struct gmx_file_position_t;
54 struct ObservablesHistory;
55 struct t_commrec;
56 struct t_fileio;
57 struct t_inputrec;
58 class t_state;
59 struct t_trxframe;
60 enum class IntegrationAlgorithm : int;
61 enum class SwapType : int;
62 enum class LambdaWeightCalculation : int;
63 enum class CheckPointVersion : int;
64
65 namespace gmx
66 {
67
68 struct MDModulesNotifiers;
69 class KeyValueTreeObject;
70 class ReadCheckpointDataHolder;
71 class WriteCheckpointDataHolder;
72
73 /*! \brief Read to a key-value-tree value used for checkpointing.
74  *
75  * \tparam ValueType
76  *
77  * \param[in] value the value to be checkpointed
78  * \param[in] name name of the value to be checkpointed
79  * \param[in] identifier uniquely identifies the module that is checkpointing
80  *                       typically the module name
81  * \param[in] kvt the key value tree to read from
82  *
83  * \throws InternalError if kvt does not contain requested value.
84  * \note Triggers assertion if value type is not correct.
85  */
86 template<typename ValueType>
87 void readKvtCheckpointValue(compat::not_null<ValueType*> value,
88                             const std::string&           name,
89                             const std::string&           identifier,
90                             const KeyValueTreeObject&    kvt);
91 //! \copydoc readKvtCheckpointValue
92 extern template void readKvtCheckpointValue(compat::not_null<std::int64_t*> value,
93                                             const std::string&              name,
94                                             const std::string&              identifier,
95                                             const KeyValueTreeObject&       kvt);
96 //! \copydoc readKvtCheckpointValue
97 extern template void readKvtCheckpointValue(compat::not_null<real*>   value,
98                                             const std::string&        name,
99                                             const std::string&        identifier,
100                                             const KeyValueTreeObject& kvt);
101
102 /*! \brief Write to a key-value-tree used for checkpointing.
103  *
104  * \tparam ValueType
105  *
106  * \param[in] value name of the value to be checkpointed
107  * \param[in] name the value to be checkpointed
108  * \param[in] identifier uniquely identifies the module that is checkpointing
109  *                       typically the module name
110  * \param[in] kvtBuilder the key-value-tree builder used to store the checkpoint values
111  */
112 template<typename ValueType>
113 void writeKvtCheckpointValue(const ValueType&          value,
114                              const std::string&        name,
115                              const std::string&        identifier,
116                              KeyValueTreeObjectBuilder kvtBuilder);
117 //! \copydoc writeKvtCheckpointValue
118 extern template void writeKvtCheckpointValue(const std::int64_t&       value,
119                                              const std::string&        name,
120                                              const std::string&        identifier,
121                                              KeyValueTreeObjectBuilder kvtBuilder);
122 //! \copydoc writeKvtCheckpointValue
123 extern template void writeKvtCheckpointValue(const real&               value,
124                                              const std::string&        name,
125                                              const std::string&        identifier,
126                                              KeyValueTreeObjectBuilder kvtBuilder);
127
128 /*! \libinternal
129  * \brief Provides the MDModules with the checkpointed data on the master rank.
130  */
131 struct MDModulesCheckpointReadingDataOnMaster
132 {
133     //! The data of the MDModules that is stored in the checkpoint file
134     const KeyValueTreeObject& checkpointedData_;
135     //! The version of the read ceckpoint file
136     CheckPointVersion checkpointFileVersion_;
137 };
138
139 /*! \libinternal
140  * \brief Provides the MDModules with the communication record to broadcast.
141  */
142 struct MDModulesCheckpointReadingBroadcast
143 {
144     //! The communicator
145     MPI_Comm communicator_;
146     //! Whether the run is executed in parallel
147     bool isParallelRun_;
148     //! The version of the read file version
149     CheckPointVersion checkpointFileVersion_;
150 };
151
152 /*! \libinternal \brief Writing the MDModules data to a checkpoint file.
153  */
154 struct MDModulesWriteCheckpointData
155 {
156     //! Builder for the Key-Value-Tree to store the MDModule checkpoint data
157     KeyValueTreeObjectBuilder builder_;
158     //! The version of the read file version
159     CheckPointVersion checkpointFileVersion_;
160 };
161
162 } // namespace gmx
163
164 /* the name of the environment variable to disable fsync failure checks with */
165 #define GMX_IGNORE_FSYNC_FAILURE_ENV "GMX_IGNORE_FSYNC_FAILURE"
166
167 // TODO Replace this mechanism with std::array<char, 1024> or similar.
168 #define CPTSTRLEN 1024
169
170 /*! \brief Enum of values that describe the contents of a cpt file
171  * whose format matches a version number
172  *
173  * The enum helps the code be more self-documenting and ensure merges
174  * do not silently resolve when two patches make the same bump. When
175  * adding new functionality, add a new element just above Count
176  * in this enumeration, and write code that does the right thing
177  * according to the value of file_version.
178  */
179 enum class CheckPointVersion : int
180 {
181     //! Unknown initial version
182     UnknownVersion0,
183     //! Unknown version 1
184     UnknownVersion1,
185     //! Store magic number to validate file.
186     AddMagicNumber,
187     //! Store which sim this is.
188     SafeSimulationPart,
189     //! Store energy data.
190     EkinDataAndFlags,
191     //! Store current step.
192     SafeSteps,
193     //! Cutoff before version 4.5.
194     Version45,
195     //! Unknown version 7.
196     UnknownVersion7,
197     //! Store checksum.
198     FileChecksumAndSize,
199     //! Unknown version 9.
200     UnknownVersion9,
201     //! Safe NH chains for thermostat.
202     NoseHooverThermostat,
203     //! Safe NH chains for barostat.
204     NoseHooverBarostat,
205     //! Safe host info.
206     HostInformation,
207     //! Activate double build.
208     DoublePrecisionBuild,
209     //! Lambda stuff.
210     LambdaStateAndHistory,
211     //! ED stuff.
212     EssentialDynamics,
213     //! Swap state stuff.
214     SwapState,
215     //! AWH history flags added.
216     AwhHistoryFlags,
217     //! Remove functionality that makes mdrun builds non-reproducible.
218     RemoveBuildMachineInformation,
219     //! Allow using COM of previous step as pull group PBC reference.
220     ComPrevStepAsPullGroupReference,
221     //! Added possibility to output average pull force and position.
222     PullAverage,
223     //! Added checkpointing for MDModules.
224     MDModules,
225     //! Added checkpointing for modular simulator.
226     ModularSimulator,
227     //! The total number of checkpoint versions.
228     Count,
229     //! Current version
230     CurrentVersion = Count - 1
231 };
232
233
234 /*!
235  * \brief
236  * Header explaining the context of a checkpoint file.
237  *
238  * TODO Expand this into being a container of all data for
239  * serialization of a checkpoint, which can be stored by the caller
240  * (e.g. so that mdrun doesn't have to open the checkpoint twice).
241  * This will separate issues of allocation from those of
242  * serialization, help separate comparison from reading, and have
243  * better defined transformation functions to/from trajectory frame
244  * data structures.
245  *
246  * Several fields were once written to checkpoint file headers, but
247  * have been removed. So that old files can continue to be read,
248  * the names of such fields contain the string "_UNUSED" so that it
249  * is clear they should not be used.
250  */
251 struct CheckpointHeaderContents
252 {
253     //! Version of checkpoint file read from disk.
254     CheckPointVersion file_version;
255     //! Version string.
256     char version[CPTSTRLEN];
257     //! Deprecated string for time.
258     char btime_UNUSED[CPTSTRLEN];
259     //! Deprecated string for user.
260     char buser_UNUSED[CPTSTRLEN];
261     //! Deprecated string for host.
262     char bhost_UNUSED[CPTSTRLEN];
263     //! Value for precision.
264     int double_prec;
265     //! Program string.
266     char fprog[CPTSTRLEN];
267     //! Time string.
268     char ftime[CPTSTRLEN];
269     //! Which integrator is in use.
270     IntegrationAlgorithm eIntegrator;
271     //! Which part of the simulation this is.
272     int simulation_part;
273     //! Which step the checkpoint is at.
274     int64_t step;
275     //! Current simulation time.
276     double t;
277     //! Number of nodes used for simulation,
278     int nnodes;
279     //! Domain decomposition settings?
280     ivec dd_nc;
281     //! Number of separate PME ranks.
282     int npme;
283     //! Number of atoms.
284     int natoms;
285     //! Number of temperature coupling groups.
286     int ngtc;
287     //! Number of Nose-Hoover pressure coupling chains.
288     int nnhpres;
289     //! Length of Nose-Hoover chains.
290     int nhchainlength;
291     //! Current FEP lambda state.
292     int nlambda;
293     //! Current state flags.
294     int flags_state;
295     //! Flags for kinetic energy.
296     int flags_eks;
297     //! Flags for energy history.
298     int flags_enh;
299     //! Flags for pull history.
300     int flagsPullHistory;
301     //! Flags for mystery history.
302     int flags_dfh;
303     //! Flags for AWH history.
304     int flags_awhh;
305     //! Essential dynamics states.
306     int nED;
307     //! Enum for coordinate swapping.
308     SwapType eSwapCoords;
309     //! Whether the checkpoint was written by modular simulator.
310     bool isModularSimulatorCheckpoint = false;
311 };
312
313 /*! \brief Low-level checkpoint writing function */
314 void write_checkpoint_data(t_fileio*                         fp,
315                            CheckpointHeaderContents          headerContents,
316                            gmx_bool                          bExpanded,
317                            LambdaWeightCalculation           elamstats,
318                            t_state*                          state,
319                            ObservablesHistory*               observablesHistory,
320                            const gmx::MDModulesNotifiers&    mdModulesNotifiers,
321                            std::vector<gmx_file_position_t>* outputfiles,
322                            gmx::WriteCheckpointDataHolder*   modularSimulatorCheckpointData);
323
324 /* Loads a checkpoint from fn for run continuation.
325  * Generates a fatal error on system size mismatch.
326  * The master node reads the file
327  * and communicates all the modified number of steps,
328  * but not the state itself.
329  * With reproducibilityRequested warns about version, build, #ranks differences.
330  */
331 void load_checkpoint(const char*                    fn,
332                      t_fileio*                      logfio,
333                      const t_commrec*               cr,
334                      const ivec                     dd_nc,
335                      t_inputrec*                    ir,
336                      t_state*                       state,
337                      ObservablesHistory*            observablesHistory,
338                      gmx_bool                       reproducibilityRequested,
339                      const gmx::MDModulesNotifiers& mdModulesNotifiers,
340                      gmx::ReadCheckpointDataHolder* modularSimulatorCheckpointData,
341                      bool                           useModularSimulator);
342
343 /* Read everything that can be stored in t_trxframe from a checkpoint file */
344 void read_checkpoint_trxframe(struct t_fileio* fp, t_trxframe* fr);
345
346 /* Print the complete contents of checkpoint file fn to out */
347 void list_checkpoint(const char* fn, FILE* out);
348
349 /*!\brief Read simulation step and part from a checkpoint file
350  *
351  * Used by tune_pme to handle tuning with a checkpoint file as part of the input.
352  *
353  * \param[in]  filename         Name of checkpoint file
354  * \param[out] simulation_part  The part of the simulation that wrote the checkpoint
355  * \param[out] step             The final step number of the simulation that wrote the checkpoint
356  *
357  * The output variables will both contain 0 if filename is NULL, the file
358  * does not exist, or is not readable. */
359 void read_checkpoint_part_and_step(const char* filename, int* simulation_part, int64_t* step);
360
361 /*!\brief Return header information from an open checkpoint file.
362  *
363  * Used by mdrun to handle restarts
364  *
365  * \param[in]  fp               Handle to open checkpoint file
366  * \param[out] outputfiles      Container of output file names from the previous run. */
367 CheckpointHeaderContents
368 read_checkpoint_simulation_part_and_filenames(t_fileio* fp, std::vector<gmx_file_position_t>* outputfiles);
369
370 #endif