c15e7321af273a20dd0f7aa3d8f79f28e858201f
[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,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 #ifndef _checkpoint_h
39 #define _checkpoint_h
40
41 #include <cstdio>
42
43 #include <vector>
44
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/keyvaluetreebuilder.h"
48
49 class energyhistory_t;
50 struct gmx_file_position_t;
51 struct ObservablesHistory;
52 struct t_commrec;
53 struct t_fileio;
54 struct t_inputrec;
55 class t_state;
56 struct t_trxframe;
57
58 namespace gmx
59 {
60
61 struct MdModulesNotifier;
62 class KeyValueTreeObject;
63
64 struct MdModulesCheckpointReadingDataOnMaster
65 {
66     //! The data of the MdModules that is stored in the checkpoint file
67     const KeyValueTreeObject& checkpointedData_;
68     //! The version of the read ceckpoint file
69     int checkpointFileVersion_;
70 };
71
72 /*! \libinternal
73  * \brief Provides the MdModules with the communication record to broadcast.
74  */
75 struct MdModulesCheckpointReadingBroadcast
76 {
77     //! The communication record
78     const t_commrec& cr_;
79     //! The version of the read file version
80     int checkpointFileVersion_;
81 };
82
83 /*! \libinternal \brief Writing the MdModules data to a checkpoint file.
84  */
85 struct MdModulesWriteCheckpointData
86 {
87     //! Builder for the Key-Value-Tree to store the MdModule checkpoint data
88     KeyValueTreeObjectBuilder builder_;
89     //! The version of the read file version
90     int checkpointFileVersion_;
91 };
92
93 } // namespace gmx
94
95 /* the name of the environment variable to disable fsync failure checks with */
96 #define GMX_IGNORE_FSYNC_FAILURE_ENV "GMX_IGNORE_FSYNC_FAILURE"
97
98 // TODO Replace this mechanism with std::array<char, 1024> or similar.
99 #define CPTSTRLEN 1024
100
101 /*!
102  * \brief
103  * Header explaining the context of a checkpoint file.
104  *
105  * TODO Expand this into being a container of all data for
106  * serialization of a checkpoint, which can be stored by the caller
107  * (e.g. so that mdrun doesn't have to open the checkpoint twice).
108  * This will separate issues of allocation from those of
109  * serialization, help separate comparison from reading, and have
110  * better defined transformation functions to/from trajectory frame
111  * data structures.
112  *
113  * Several fields were once written to checkpoint file headers, but
114  * have been removed. So that old files can continue to be read,
115  * the names of such fields contain the string "_UNUSED" so that it
116  * is clear they should not be used.
117  */
118 struct CheckpointHeaderContents
119 {
120     //! Version of checkpoint file read from disk.
121     int file_version;
122     //! Version string.
123     char version[CPTSTRLEN];
124     //! Deprecated string for time.
125     char btime_UNUSED[CPTSTRLEN];
126     //! Deprecated string for user.
127     char buser_UNUSED[CPTSTRLEN];
128     //! Deprecated string for host.
129     char bhost_UNUSED[CPTSTRLEN];
130     //! Value for precision.
131     int double_prec;
132     //! Program string.
133     char fprog[CPTSTRLEN];
134     //! Time string.
135     char ftime[CPTSTRLEN];
136     //! Which integrator is in use.
137     int eIntegrator;
138     //! Which part of the simulation this is.
139     int simulation_part;
140     //! Which step the checkpoint is at.
141     int64_t step;
142     //! Current simulation time.
143     double t;
144     //! Number of nodes used for simulation,
145     int nnodes;
146     //! Domain decomposition settings?
147     ivec dd_nc;
148     //! Number of separate PME ranks.
149     int npme;
150     //! Number of atoms.
151     int natoms;
152     //! Number of temperature coupling groups.
153     int ngtc;
154     //! Number of Nose-Hoover pressure coupling chains.
155     int nnhpres;
156     //! Length of Nose-Hoover chains.
157     int nhchainlength;
158     //! Current FEP lambda state.
159     int nlambda;
160     //! Current state flags.
161     int flags_state;
162     //! Flags for kinetic energy.
163     int flags_eks;
164     //! Flags for energy history.
165     int flags_enh;
166     //! Flags for pull history.
167     int flagsPullHistory;
168     //! Flags for mystery history.
169     int flags_dfh;
170     //! Flags for AWH history.
171     int flags_awhh;
172     //! Essential dynamics states.
173     int nED;
174     //! Enum for coordinate swapping.
175     int eSwapCoords;
176 };
177
178 /* Write a checkpoint to <fn>.cpt
179  * Appends the _step<step>.cpt with bNumberAndKeep,
180  * otherwise moves the previous <fn>.cpt to <fn>_prev.cpt
181  */
182 void write_checkpoint(const char*                   fn,
183                       gmx_bool                      bNumberAndKeep,
184                       FILE*                         fplog,
185                       const t_commrec*              cr,
186                       ivec                          domdecCells,
187                       int                           nppnodes,
188                       int                           eIntegrator,
189                       int                           simulation_part,
190                       gmx_bool                      bExpanded,
191                       int                           elamstats,
192                       int64_t                       step,
193                       double                        t,
194                       t_state*                      state,
195                       ObservablesHistory*           observablesHistory,
196                       const gmx::MdModulesNotifier& notifier);
197
198 /* Loads a checkpoint from fn for run continuation.
199  * Generates a fatal error on system size mismatch.
200  * The master node reads the file
201  * and communicates all the modified number of steps,
202  * but not the state itself.
203  * With reproducibilityRequested warns about version, build, #ranks differences.
204  */
205 void load_checkpoint(const char*                   fn,
206                      t_fileio*                     logfio,
207                      const t_commrec*              cr,
208                      const ivec                    dd_nc,
209                      t_inputrec*                   ir,
210                      t_state*                      state,
211                      ObservablesHistory*           observablesHistory,
212                      gmx_bool                      reproducibilityRequested,
213                      const gmx::MdModulesNotifier& mdModulesNotifier);
214
215 /* Read everything that can be stored in t_trxframe from a checkpoint file */
216 void read_checkpoint_trxframe(struct t_fileio* fp, t_trxframe* fr);
217
218 /* Print the complete contents of checkpoint file fn to out */
219 void list_checkpoint(const char* fn, FILE* out);
220
221 /*!\brief Read simulation step and part from a checkpoint file
222  *
223  * Used by tune_pme to handle tuning with a checkpoint file as part of the input.
224  *
225  * \param[in]  filename         Name of checkpoint file
226  * \param[out] simulation_part  The part of the simulation that wrote the checkpoint
227  * \param[out] step             The final step number of the simulation that wrote the checkpoint
228  *
229  * The output variables will both contain 0 if filename is NULL, the file
230  * does not exist, or is not readable. */
231 void read_checkpoint_part_and_step(const char* filename, int* simulation_part, int64_t* step);
232
233 /*!\brief Return header information from an open checkpoint file.
234  *
235  * Used by mdrun to handle restarts
236  *
237  * \param[in]  fp               Handle to open checkpoint file
238  * \param[out] outputfiles      Container of output file names from the previous run. */
239 CheckpointHeaderContents
240 read_checkpoint_simulation_part_and_filenames(t_fileio* fp, std::vector<gmx_file_position_t>* outputfiles);
241
242 #endif