Warn for type mismatch for gmx printf like functions 2/3
[gromacs.git] / src / gromacs / utility / keyvaluetree.h
blob293289d663ab49777fff640f98a546b457f269f0
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018, 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
37 * Declares a data structure for JSON-like structured key-value mapping.
39 * A tree is composed of nodes that can have different types:
40 * - _Value_ (gmx::KeyValueTreeValue) is a generic node that can
41 * represent either a scalar value of arbitrary type, or an object or
42 * an array.
43 * - _Array_ (gmx::KeyValueTreeArray) is a collection of any number of values
44 * (including zero). The values can be of any type and different types
45 * can be mixed in the same array.
46 * - _Object_ (gmx::KeyValueTreeObject) is a collection of properties.
47 * Each property must have a unique key. Order of properties is preserved,
48 * i.e., they can be iterated in the order they were added.
49 * - _Property_ (gmx::KeyValueTreeProperty) is an arbitrary type of value
50 * associated with a string key.
51 * The root object of a tree is typically an object, but could also be an
52 * array. The data structure itself does not enforce any other constraints,
53 * but the context in which it is used can limit the allowed scalar types or,
54 * e.g., require arrays to have values of uniform type. Also, several
55 * operations defined for the structure (string output, comparison,
56 * serialization, etc.) only work on a limited set of scalar types, or have
57 * limitations with the types of trees they work on (in particular, arrays are
58 * currently poorly supported).
60 * \author Teemu Murtola <teemu.murtola@gmail.com>
61 * \inlibraryapi
62 * \ingroup module_utility
64 #ifndef GMX_UTILITY_KEYVALUETREE_H
65 #define GMX_UTILITY_KEYVALUETREE_H
67 #include <algorithm>
68 #include <functional>
69 #include <map>
70 #include <string>
71 #include <utility>
72 #include <vector>
74 #include "gromacs/utility/real.h"
75 #include "gromacs/utility/variant.h"
77 namespace gmx
80 class KeyValueTreeArray;
81 class KeyValueTreeObject;
82 class TextWriter;
84 /*! \libinternal \brief
85 * Identifies an entry in a key-value tree.
87 * This class is mainly an internal utility within the key-value tree
88 * implementation, but it is exposed on the API level where string-based
89 * specification of a location in the tree is necessary. Constructors are not
90 * explicit to allow passing a simple string in contexts where a
91 * KeyValueTreePath is expected.
93 * The string specifying a location should start with a `/`, followed by the
94 * names of the properties separated by `/`. For example, `/a/b/c` specifies
95 * property `c` in an object that is the value of `b` in an object that is the
96 * value of `a` at the root of the tree.
97 * Currently, there is no support for specifying paths to values within arrays
98 * (since none of the places where this is used implement array handling,
99 * either).
101 * \inlibraryapi
102 * \ingroup module_utility
104 class KeyValueTreePath
106 public:
107 //! Creates an empty path (corresponds to the root object).
108 KeyValueTreePath() = default;
109 //! Creates a path from given string representation.
110 KeyValueTreePath(const char *path);
111 //! Creates a path from given string representation.
112 KeyValueTreePath(const std::string &path);
114 //! Adds another element to the path, making it a child of the old path.
115 void append(const std::string &key) { path_.push_back(key); }
116 //! Adds elements from another path to the path.
117 void append(const KeyValueTreePath &other)
119 auto elements = other.elements();
120 path_.insert(path_.end(), elements.begin(), elements.end());
122 //! Removes the last element in the path, making it the parent path.
123 void pop_back() { return path_.pop_back(); }
124 //! Removes and returns the last element in the path.
125 std::string pop_last()
127 std::string result = std::move(path_.back());
128 path_.pop_back();
129 return result;
132 //! Whether the path is empty (pointing to the root object).
133 bool empty() const { return path_.empty(); }
134 //! Returns the number of elements (=nesting level) in the path.
135 size_t size() const { return path_.size(); }
136 //! Returns the i'th path element.
137 const std::string &operator[](int i) const { return path_[i]; }
138 //! Returns all the path elements.
139 const std::vector<std::string> &elements() const { return path_; }
141 //! Formats the path as a string for display.
142 std::string toString() const;
144 private:
145 std::vector<std::string> path_;
148 //! \cond libapi
150 //! Combines two paths as with KeyValueTreePath::append().
151 inline KeyValueTreePath operator+(const KeyValueTreePath &a, const KeyValueTreePath &b)
153 KeyValueTreePath result(a);
154 result.append(b);
155 return result;
158 //! Combines an element to a path as with KeyValueTreePath::append().
159 inline KeyValueTreePath operator+(const KeyValueTreePath &a, const std::string &b)
161 KeyValueTreePath result(a);
162 result.append(b);
163 return result;
165 //! \endcond
167 class KeyValueTreeValue
169 public:
170 //! Returns whether the value is an array (KeyValueTreeArray).
171 bool isArray() const;
172 //! Returns whether the value is an object (KeyValueTreeObject).
173 bool isObject() const;
174 //! Returns whether the value is of a given type.
175 template <typename T>
176 bool isType() const { return value_.isType<T>(); }
177 //! Returns the type of the value.
178 std::type_index type() const { return value_.type(); }
180 KeyValueTreeArray &asArray();
181 KeyValueTreeObject &asObject();
182 const KeyValueTreeArray &asArray() const;
183 const KeyValueTreeObject &asObject() const;
184 template <typename T>
185 const T &cast() const { return value_.cast<T>(); }
187 //! Returns the raw Variant value (always possible).
188 const Variant &asVariant() const { return value_; }
190 private:
191 explicit KeyValueTreeValue(Variant &&value) : value_(std::move(value)) {}
193 Variant value_;
195 friend class KeyValueTreeBuilder;
196 friend class KeyValueTreeObjectBuilder;
197 friend class KeyValueTreeValueBuilder;
200 class KeyValueTreeArray
202 public:
203 //! Whether all elements of the array are objects.
204 bool isObjectArray() const
206 return std::all_of(values_.begin(), values_.end(),
207 std::mem_fn(&KeyValueTreeValue::isObject));
210 //! Returns the values in the array.
211 const std::vector<KeyValueTreeValue> &values() const { return values_; }
213 private:
214 std::vector<KeyValueTreeValue> values_;
216 friend class KeyValueTreeArrayBuilderBase;
219 class KeyValueTreeProperty
221 public:
222 const std::string &key() const { return value_->first; }
223 const KeyValueTreeValue &value() const { return value_->second; }
225 private:
226 typedef std::map<std::string, KeyValueTreeValue>::const_iterator
227 IteratorType;
229 explicit KeyValueTreeProperty(IteratorType value) : value_(value) {}
231 IteratorType value_;
233 friend class KeyValueTreeObject;
234 friend class KeyValueTreeObjectBuilder;
237 class KeyValueTreeObject
239 public:
240 KeyValueTreeObject() = default;
241 //! Creates a deep copy of an object.
242 KeyValueTreeObject(const KeyValueTreeObject &other)
244 for (const auto &value : other.values_)
246 auto iter = valueMap_.insert(std::make_pair(value.key(), value.value())).first;
247 values_.push_back(KeyValueTreeProperty(iter));
250 //! Assigns a deep copy of an object.
251 KeyValueTreeObject &operator=(const KeyValueTreeObject &other)
253 KeyValueTreeObject tmp(other);
254 std::swap(tmp.valueMap_, valueMap_);
255 std::swap(tmp.values_, values_);
256 return *this;
258 //! Default move constructor.
259 //NOLINTNEXTLINE(performance-noexcept-move-constructor) bug #38733
260 KeyValueTreeObject(KeyValueTreeObject &&) = default;
261 //! Default move assignment.
262 KeyValueTreeObject &operator=(KeyValueTreeObject &&) = default;
264 /*! \brief
265 * Returns all properties in the object.
267 * The properties are in the order they were added to the object.
269 const std::vector<KeyValueTreeProperty> &properties() const { return values_; }
271 //! Whether a property with given key exists.
272 bool keyExists(const std::string &key) const
274 return valueMap_.find(key) != valueMap_.end();
276 //! Returns value for a given key.
277 const KeyValueTreeValue &operator[](const std::string &key) const
279 GMX_ASSERT(keyExists(key), "Accessing non-existent value");
280 return valueMap_.at(key);
283 /*! \brief
284 * Returns whether the given object shares any keys with `this`.
286 bool hasDistinctProperties(const KeyValueTreeObject &obj) const;
288 private:
289 //! Keeps the properties by key.
290 std::map<std::string, KeyValueTreeValue> valueMap_;
291 //! Keeps the insertion order of properties.
292 std::vector<KeyValueTreeProperty> values_;
294 friend class KeyValueTreeObjectBuilder;
297 /********************************************************************
298 * Inline functions that could not be declared within the classes
301 inline bool KeyValueTreeValue::isArray() const
303 return value_.isType<KeyValueTreeArray>();
305 inline bool KeyValueTreeValue::isObject() const
307 return value_.isType<KeyValueTreeObject>();
309 inline const KeyValueTreeArray &KeyValueTreeValue::asArray() const
311 return value_.cast<KeyValueTreeArray>();
313 inline const KeyValueTreeObject &KeyValueTreeValue::asObject() const
315 return value_.cast<KeyValueTreeObject>();
317 inline KeyValueTreeArray &KeyValueTreeValue::asArray()
319 return value_.castRef<KeyValueTreeArray>();
321 inline KeyValueTreeObject &KeyValueTreeValue::asObject()
323 return value_.castRef<KeyValueTreeObject>();
326 //! \cond libapi
327 /*! \brief
328 * Writes a human-readable representation of the tree with given writer.
330 * The output format is designed to be readable by humans; if some
331 * particular machine-readable format is needed, that should be
332 * implemented outside the generic key-value tree code.
334 * \ingroup module_utility
336 void dumpKeyValueTree(TextWriter *writer, const KeyValueTreeObject &tree);
338 /*! \brief
339 * Compares two KeyValueTrees and prints any differences.
341 * \ingroup module_utility
343 void compareKeyValueTrees(TextWriter *writer,
344 const KeyValueTreeObject &tree1,
345 const KeyValueTreeObject &tree2,
346 real ftol,
347 real abstol);
349 //! Helper function to format a simple KeyValueTreeValue.
350 static inline std::string
351 simpleValueToString(const KeyValueTreeValue &value)
353 return simpleValueToString(value.asVariant());
356 //! \endcond
358 } // namespace gmx
360 #endif