Merge branch release-2021
[alexxy/gromacs.git] / src / gromacs / mdtypes / checkpointdata.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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.
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 Provides the checkpoint data structure for the modular simulator
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_mdtypes
40  */
41
42 #ifndef GMX_MODULARSIMULATOR_CHECKPOINTDATA_H
43 #define GMX_MODULARSIMULATOR_CHECKPOINTDATA_H
44
45 #include <optional>
46
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/utility/arrayref.h"
49 #include "gromacs/utility/exceptions.h"
50 #include "gromacs/utility/keyvaluetreebuilder.h"
51 #include "gromacs/utility/stringutil.h"
52
53 namespace gmx
54 {
55 class ISerializer;
56
57 /*! \libinternal
58  * \brief The operations on CheckpointData
59  *
60  * This enum defines the two modes of operation on CheckpointData objects,
61  * reading and writing. This allows to template all access functions, which
62  * in turn enables clients to write a single function for read and write
63  * access, eliminating the risk of having read and write functions getting
64  * out of sync.
65  *
66  * \ingroup module_modularsimulator
67  */
68 enum class CheckpointDataOperation
69 {
70     Read,
71     Write,
72     Count
73 };
74
75 /*! \internal
76  * \brief Get an ArrayRef whose const-ness is defined by the checkpointing operation
77  *
78  * \tparam operation  Whether we are reading or writing
79  * \tparam T          The type of values stored in the ArrayRef
80  * \param container   The container the ArrayRef is referencing to
81  * \return            The ArrayRef
82  *
83  * \see ArrayRef
84  *
85  * \ingroup module_modularsimulator
86  */
87 template<CheckpointDataOperation operation, typename T>
88 ArrayRef<std::conditional_t<operation == CheckpointDataOperation::Write || std::is_const<T>::value, const typename T::value_type, typename T::value_type>>
89 makeCheckpointArrayRef(T& container)
90 {
91     return container;
92 }
93
94 /*! \internal
95  * \ingroup module_modularsimulator
96  * \brief Struct allowing to check if data is serializable through the KeyValueTree serializer
97  *
98  * This list of types is copied from ValueSerializer::initSerializers()
99  * Having this here allows us to catch errors at compile time
100  * instead of having cryptic runtime errors
101  */
102 template<typename T>
103 struct IsSerializableType
104 {
105     static bool const value = std::is_same<T, std::string>::value || std::is_same<T, bool>::value
106                               || std::is_same<T, int>::value || std::is_same<T, int64_t>::value
107                               || std::is_same<T, float>::value || std::is_same<T, double>::value;
108 };
109
110 /*! \internal
111  * \ingroup module_modularsimulator
112  * \brief Struct allowing to check if enum has a serializable underlying type
113  */
114 //! {
115 template<typename T, bool = std::is_enum<T>::value>
116 struct IsSerializableEnum
117 {
118     static bool const value = IsSerializableType<std::underlying_type_t<T>>::value;
119 };
120 template<typename T>
121 struct IsSerializableEnum<T, false>
122 {
123     static bool const value = false;
124 };
125 //! }
126
127 /*! \libinternal
128  * \ingroup module_modularsimulator
129  * \brief Data type hiding checkpoint implementation details
130  *
131  * This data type allows to separate the implementation details of the
132  * checkpoint writing / reading from the implementation of the checkpoint
133  * clients. Checkpoint clients interface via the methods of the CheckpointData
134  * object, and do not need knowledge of data types used to store the data.
135  *
136  * Templating allows checkpoint clients to have symmetric (templated)
137  * implementations for checkpoint reading and writing.
138  *
139  * CheckpointData objects are dispatched via [Write|Read]CheckpointDataHolder
140  * objects, which interact with the checkpoint reading from / writing to
141  * file.
142  */
143
144 template<CheckpointDataOperation operation>
145 class CheckpointData;
146
147 //! Convenience shortcut for reading checkpoint data.
148 using ReadCheckpointData = CheckpointData<CheckpointDataOperation::Read>;
149 //! Convenience shortcut for writing checkpoint data.
150 using WriteCheckpointData = CheckpointData<CheckpointDataOperation::Write>;
151
152 template<>
153 class CheckpointData<CheckpointDataOperation::Read>
154 {
155 public:
156     /*! \brief Read or write a single value from / to checkpoint
157      *
158      * Allowed scalar types include std::string, bool, int, int64_t,
159      * float, double, or any enum with one of the previously mentioned
160      * scalar types as underlying type. Type compatibility is checked
161      * at compile time.
162      *
163      * \tparam operation  Whether we are reading or writing
164      * \tparam T          The type of the value
165      * \param key         The key to [read|write] the value [from|to]
166      * \param value       The value to [read|write]
167      */
168     //! {
169     template<typename T>
170     std::enable_if_t<IsSerializableType<T>::value, void> scalar(const std::string& key, T* value) const;
171     template<typename T>
172     std::enable_if_t<IsSerializableEnum<T>::value, void> enumScalar(const std::string& key, T* value) const;
173     //! }
174
175     /*! \brief Read or write an ArrayRef from / to checkpoint
176      *
177      * Allowed types stored in the ArrayRef include std::string, bool, int,
178      * int64_t, float, double, and gmx::RVec. Type compatibility is checked
179      * at compile time.
180      *
181      * \tparam operation  Whether we are reading or writing
182      * \tparam T          The type of values stored in the ArrayRef
183      * \param key         The key to [read|write] the ArrayRef [from|to]
184      * \param values      The ArrayRef to [read|write]
185      */
186     //! {
187     // Read ArrayRef of scalar
188     template<typename T>
189     std::enable_if_t<IsSerializableType<T>::value, void> arrayRef(const std::string& key,
190                                                                   ArrayRef<T>        values) const;
191     // Read ArrayRef of RVec
192     void arrayRef(const std::string& key, ArrayRef<RVec> values) const;
193     //! }
194
195     /*! \brief Read or write a tensor from / to checkpoint
196      *
197      * \tparam operation  Whether we are reading or writing
198      * \param key         The key to [read|write] the tensor [from|to]
199      * \param values      The tensor to [read|write]
200      */
201     void tensor(const std::string& key, ::tensor values) const;
202
203     /*! \brief Return a subset of the current CheckpointData
204      *
205      * \tparam operation  Whether we are reading or writing
206      * \param key         The key to [read|write] the sub data [from|to]
207      * \return            A CheckpointData object representing a subset of the current object
208      */
209     //!{
210     CheckpointData subCheckpointData(const std::string& key) const;
211     //!}
212
213 private:
214     //! KV tree read from checkpoint
215     const KeyValueTreeObject* inputTree_ = nullptr;
216
217     //! Construct an input checkpoint data object
218     explicit CheckpointData(const KeyValueTreeObject& inputTree);
219
220     // Only holders should build
221     friend class ReadCheckpointDataHolder;
222 };
223
224 template<>
225 class CheckpointData<CheckpointDataOperation::Write>
226 {
227 public:
228     //! \copydoc CheckpointData<CheckpointDataOperation::Read>::scalar
229     //! {
230     template<typename T>
231     std::enable_if_t<IsSerializableType<T>::value, void> scalar(const std::string& key, const T* value);
232     template<typename T>
233     std::enable_if_t<IsSerializableEnum<T>::value, void> enumScalar(const std::string& key, const T* value);
234     //! }
235
236     //! \copydoc CheckpointData<CheckpointDataOperation::Read>::arrayRef
237     //! {
238     // Write ArrayRef of scalar
239     template<typename T>
240     std::enable_if_t<IsSerializableType<T>::value, void> arrayRef(const std::string& key,
241                                                                   ArrayRef<const T>  values);
242     // Write ArrayRef of RVec
243     void arrayRef(const std::string& key, ArrayRef<const RVec> values);
244     //! }
245
246     //! \copydoc CheckpointData<CheckpointDataOperation::Read>::tensor
247     void tensor(const std::string& key, const ::tensor values);
248
249     //! \copydoc CheckpointData<CheckpointDataOperation::Read>::subCheckpointData
250     CheckpointData subCheckpointData(const std::string& key);
251
252 private:
253     //! Builder for the tree to be written to checkpoint
254     std::optional<KeyValueTreeObjectBuilder> outputTreeBuilder_ = std::nullopt;
255
256     //! Construct an output checkpoint data object
257     explicit CheckpointData(KeyValueTreeObjectBuilder&& outputTreeBuilder);
258
259     // Only holders should build
260     friend class WriteCheckpointDataHolder;
261 };
262
263 /*! \brief Read a checkpoint version enum variable
264  *
265  * This reads the checkpoint version from file. The read version is returned.
266  *
267  * If the read version is more recent than the code version, this throws an error, since
268  * we cannot know what has changed in the meantime. Using newer checkpoint files with
269  * old code is not a functionality we can offer. Note, however, that since the checkpoint
270  * version is saved by module, older checkpoint files of all simulations that don't use
271  * that specific module can still be used.
272  *
273  * Allowing backwards compatibility of files (i.e., reading an older checkpoint file with
274  * a newer version of the code) is in the responsibility of the caller module. They can
275  * use the returned file checkpoint version to do that:
276  *
277  *     const auto fileVersion = checkpointVersion(checkpointData, "version", c_currentVersion);
278  *     if (fileVersion >= CheckpointVersion::AddedX)
279  *     {
280  *         checkpointData->scalar("x", &x_));
281  *     }
282  *
283  * @tparam VersionEnum     The type of the checkpoint version enum
284  * @param  checkpointData  A reading checkpoint data object
285  * @param  key             The key under which the version is saved - also used for error output
286  * @param  programVersion  The checkpoint version of the current code
287  * @return                 The checkpoint version read from file
288  */
289 template<typename VersionEnum>
290 VersionEnum checkpointVersion(const ReadCheckpointData* checkpointData,
291                               const std::string&        key,
292                               const VersionEnum         programVersion)
293 {
294     VersionEnum fileVersion;
295     checkpointData->enumScalar(key, &fileVersion);
296     if (fileVersion > programVersion)
297     {
298         throw FileIOError(
299                 formatString("The checkpoint file contains a %s that is more recent than the "
300                              "current program version and is not backward compatible.",
301                              key.c_str()));
302     }
303     return fileVersion;
304 }
305
306 /*! \brief Write the current code checkpoint version enum variable
307  *
308  * Write the current program checkpoint version to the checkpoint data object.
309  * Returns the written checkpoint version to mirror the signature of the reading version.
310  *
311  * @tparam VersionEnum     The type of the checkpoint version enum
312  * @param  checkpointData  A writing checkpoint data object
313  * @param  key             The key under which the version is saved
314  * @param  programVersion  The checkpoint version of the current code
315  * @return                 The checkpoint version written to file
316  */
317 template<typename VersionEnum>
318 VersionEnum checkpointVersion(WriteCheckpointData* checkpointData,
319                               const std::string&   key,
320                               const VersionEnum    programVersion)
321 {
322     checkpointData->enumScalar(key, &programVersion);
323     return programVersion;
324 }
325
326 inline ReadCheckpointData::CheckpointData(const KeyValueTreeObject& inputTree) :
327     inputTree_(&inputTree)
328 {
329 }
330
331 inline WriteCheckpointData::CheckpointData(KeyValueTreeObjectBuilder&& outputTreeBuilder) :
332     outputTreeBuilder_(outputTreeBuilder)
333 {
334 }
335 /*! \libinternal
336  * \brief Holder for read checkpoint data
337  *
338  * A ReadCheckpointDataHolder object is passed to the checkpoint reading
339  * functionality, and then passed into the SimulatorBuilder object. It
340  * holds the KV-tree read from file and dispatches CheckpointData objects
341  * to the checkpoint clients.
342  */
343 class ReadCheckpointDataHolder
344 {
345 public:
346     //! Check whether a key exists
347     [[nodiscard]] bool keyExists(const std::string& key) const;
348
349     //! Return vector of existing keys
350     [[nodiscard]] std::vector<std::string> keys() const;
351
352     //! Deserialize serializer content into the CheckpointData object
353     void deserialize(ISerializer* serializer);
354
355     /*! \brief Return a subset of the current CheckpointData
356      *
357      * \param key         The key to [read|write] the sub data [from|to]
358      * \return            A CheckpointData object representing a subset of the current object
359      */
360     [[nodiscard]] ReadCheckpointData checkpointData(const std::string& key) const;
361
362     //! Write the contents of the Checkpoint to file
363     void dump(FILE* out) const;
364
365 private:
366     //! KV-tree read from checkpoint
367     KeyValueTreeObject checkpointTree_;
368 };
369
370 /*! \libinternal
371  * \brief Holder for write checkpoint data
372  *
373  * The WriteCheckpointDataHolder object holds the KV-tree builder and
374  * dispatches CheckpointData objects to the checkpoint clients to save
375  * their respective data. It is then passed to the checkpoint writing
376  * functionality.
377  */
378 class WriteCheckpointDataHolder
379 {
380 public:
381     //! Serialize the content of the CheckpointData object
382     void serialize(ISerializer* serializer);
383
384     /*! \brief Return a subset of the current CheckpointData
385      *
386      * \param key         The key to [read|write] the sub data [from|to]
387      * \return            A CheckpointData object representing a subset of the current object
388      */
389     [[nodiscard]] WriteCheckpointData checkpointData(const std::string& key);
390
391     /*! \brief
392      */
393     [[nodiscard]] bool empty() const;
394
395 private:
396     //! KV-tree builder
397     KeyValueTreeBuilder outputTreeBuilder_;
398     //! Whether any checkpoint data object has been requested
399     bool hasCheckpointDataBeenRequested_ = false;
400 };
401
402 // Function definitions - here to avoid template-related linker problems
403 // doxygen doesn't like these...
404 //! \cond
405 template<typename T>
406 std::enable_if_t<IsSerializableType<T>::value, void> ReadCheckpointData::scalar(const std::string& key,
407                                                                                 T* value) const
408 {
409     GMX_RELEASE_ASSERT(inputTree_, "No input checkpoint data available.");
410     *value = (*inputTree_)[key].cast<T>();
411 }
412
413 template<typename T>
414 std::enable_if_t<IsSerializableEnum<T>::value, void> ReadCheckpointData::enumScalar(const std::string& key,
415                                                                                     T* value) const
416 {
417     GMX_RELEASE_ASSERT(inputTree_, "No input checkpoint data available.");
418     std::underlying_type_t<T> castValue;
419     castValue = (*inputTree_)[key].cast<std::underlying_type_t<T>>();
420     *value    = static_cast<T>(castValue);
421 }
422
423 template<typename T>
424 inline std::enable_if_t<IsSerializableType<T>::value, void>
425 WriteCheckpointData::scalar(const std::string& key, const T* value)
426 {
427     GMX_RELEASE_ASSERT(outputTreeBuilder_, "No output checkpoint data available.");
428     outputTreeBuilder_->addValue(key, *value);
429 }
430
431 template<typename T>
432 inline std::enable_if_t<IsSerializableEnum<T>::value, void>
433 WriteCheckpointData::enumScalar(const std::string& key, const T* value)
434 {
435     GMX_RELEASE_ASSERT(outputTreeBuilder_, "No output checkpoint data available.");
436     auto castValue = static_cast<std::underlying_type_t<T>>(*value);
437     outputTreeBuilder_->addValue(key, castValue);
438 }
439
440 template<typename T>
441 inline std::enable_if_t<IsSerializableType<T>::value, void>
442 ReadCheckpointData::arrayRef(const std::string& key, ArrayRef<T> values) const
443 {
444     GMX_RELEASE_ASSERT(inputTree_, "No input checkpoint data available.");
445     GMX_RELEASE_ASSERT(values.size() >= (*inputTree_)[key].asArray().values().size(),
446                        "Read vector does not fit in passed ArrayRef.");
447     auto outputIt  = values.begin();
448     auto inputIt   = (*inputTree_)[key].asArray().values().begin();
449     auto outputEnd = values.end();
450     auto inputEnd  = (*inputTree_)[key].asArray().values().end();
451     for (; outputIt != outputEnd && inputIt != inputEnd; outputIt++, inputIt++)
452     {
453         *outputIt = inputIt->cast<T>();
454     }
455 }
456
457 template<typename T>
458 inline std::enable_if_t<IsSerializableType<T>::value, void>
459 WriteCheckpointData::arrayRef(const std::string& key, ArrayRef<const T> values)
460 {
461     GMX_RELEASE_ASSERT(outputTreeBuilder_, "No output checkpoint data available.");
462     auto builder = outputTreeBuilder_->addUniformArray<T>(key);
463     for (const auto& value : values)
464     {
465         builder.addValue(value);
466     }
467 }
468
469 inline void ReadCheckpointData::arrayRef(const std::string& key, ArrayRef<RVec> values) const
470 {
471     GMX_RELEASE_ASSERT(values.size() >= (*inputTree_)[key].asArray().values().size(),
472                        "Read vector does not fit in passed ArrayRef.");
473     auto outputIt  = values.begin();
474     auto inputIt   = (*inputTree_)[key].asArray().values().begin();
475     auto outputEnd = values.end();
476     auto inputEnd  = (*inputTree_)[key].asArray().values().end();
477     for (; outputIt != outputEnd && inputIt != inputEnd; outputIt++, inputIt++)
478     {
479         auto storedRVec = inputIt->asObject()["RVec"].asArray().values();
480         *outputIt = { storedRVec[XX].cast<real>(), storedRVec[YY].cast<real>(), storedRVec[ZZ].cast<real>() };
481     }
482 }
483
484 inline void WriteCheckpointData::arrayRef(const std::string& key, ArrayRef<const RVec> values)
485 {
486     auto builder = outputTreeBuilder_->addObjectArray(key);
487     for (const auto& value : values)
488     {
489         auto subbuilder = builder.addObject();
490         subbuilder.addUniformArray("RVec", { value[XX], value[YY], value[ZZ] });
491     }
492 }
493
494 inline void ReadCheckpointData::tensor(const std::string& key, ::tensor values) const
495 {
496     auto array     = (*inputTree_)[key].asArray().values();
497     values[XX][XX] = array[0].cast<real>();
498     values[XX][YY] = array[1].cast<real>();
499     values[XX][ZZ] = array[2].cast<real>();
500     values[YY][XX] = array[3].cast<real>();
501     values[YY][YY] = array[4].cast<real>();
502     values[YY][ZZ] = array[5].cast<real>();
503     values[ZZ][XX] = array[6].cast<real>();
504     values[ZZ][YY] = array[7].cast<real>();
505     values[ZZ][ZZ] = array[8].cast<real>();
506 }
507
508 inline void WriteCheckpointData::tensor(const std::string& key, const ::tensor values)
509 {
510     auto builder = outputTreeBuilder_->addUniformArray<real>(key);
511     builder.addValue(values[XX][XX]);
512     builder.addValue(values[XX][YY]);
513     builder.addValue(values[XX][ZZ]);
514     builder.addValue(values[YY][XX]);
515     builder.addValue(values[YY][YY]);
516     builder.addValue(values[YY][ZZ]);
517     builder.addValue(values[ZZ][XX]);
518     builder.addValue(values[ZZ][YY]);
519     builder.addValue(values[ZZ][ZZ]);
520 }
521
522 inline ReadCheckpointData ReadCheckpointData::subCheckpointData(const std::string& key) const
523 {
524     return CheckpointData((*inputTree_)[key].asObject());
525 }
526
527 inline WriteCheckpointData WriteCheckpointData::subCheckpointData(const std::string& key)
528 {
529     return CheckpointData(outputTreeBuilder_->addObject(key));
530 }
531
532 //! \endcond
533
534 } // namespace gmx
535
536 #endif // GMX_MODULARSIMULATOR_CHECKPOINTDATA_H