Add helper function to make a vector from an ArrayRef
[gromacs.git] / src / gromacs / utility / arrayref.h
blobbc1191da16fc0f2dbdb43a532a122d172ec48959
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,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 /*! \file
36 * \brief
37 * Declares gmx::ArrayRef
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \inpublicapi
41 * \ingroup module_utility
43 #ifndef GMX_UTILITY_ARRAYREF_H
44 #define GMX_UTILITY_ARRAYREF_H
46 #include <cstddef>
48 #include <array>
49 #include <iterator>
50 #include <stdexcept>
51 #include <utility>
52 #include <vector>
54 #include "gromacs/utility/gmxassert.h"
56 namespace gmx
59 /*! \brief
60 * Tag type to initialize empty array references.
62 * This type (together with appropriate constructors in ArrayRef)
63 * allows initializing any array reference to an empty value
64 * without explicitly specifying its type. This is convenient when calling
65 * a function that takes an array reference, where constructing an empty
66 * reference explicitly would otherwise require specifying the full array
67 * reference type, including the template parameter.
69 struct EmptyArrayRef {};
71 /*! \brief STL-like container for an interface to a C array of T (or part
72 * of a std::vector<T, A> or std::array<T>).
74 * \tparam T Value type of elements.
76 * This class provides an interface similar to \c std::vector<T, A>, with the
77 * following main differences:
78 * - This class does not have its own storage. Instead, it references an
79 * existing array of values (either a C-style array or part of an existing
80 * std::vector<T, A> or std::array<T>).
81 * - It is only possible to modify the values themselves through ArrayRef;
82 * it is not possible to add or remove values.
83 * - Copying objects of this type is cheap, and the copies behave identically
84 * to the original object: the copy references the same set of values.
86 * This class is useful for writing wrappers that expose a view of the
87 * internal data stored as a single vector/array, which can be a whole
88 * or part of the underlying storage.
90 * Methods in this class do not throw, except where indicated.
92 * Note that due to a Doxygen limitation, the constructor that takes a C array
93 * whose size is known at compile time does not appear in the documentation.
95 * To refer to const data of type T, ArrayRef<const T> is used. For both const
96 * and non-const std::vector and std::array an ArrayRef view can be created.
97 * Attempting to create a non-const ArrayRef of a const vector/array will result
98 * in a compiler error in the respective constructor.
100 * For SIMD types there is template specialization available
101 * (e.g. ArrayRef<SimdReal>) in gromacs/simd/simd_memory.h which should have
102 * the same functionality as much as possible.
104 * \todo
105 * This class is not complete. There are likely also methods missing (not
106 * required for current usage).
108 * \inpublicapi
109 * \ingroup module_utility
111 template <typename T>
112 class ArrayRef
114 public:
115 //! Type of values stored in the container.
116 typedef T value_type;
117 //! Type for representing size of the container.
118 typedef size_t size_type;
119 //! Type for representing difference between two container indices.
120 typedef ptrdiff_t difference_type;
121 //! Const reference to a container element.
122 typedef const T &const_reference;
123 //! Const pointer to a container element.
124 typedef const T *const_pointer;
125 //! Const iterator type for the container.
126 typedef const T *const_iterator;
127 //! Reference to a container element.
128 typedef T &reference;
129 //! Pointer to a container element.
130 typedef T *pointer;
131 //! Iterator type for the container.
132 typedef T *iterator;
133 //! Standard reverse iterator.
134 typedef std::reverse_iterator<iterator> reverse_iterator;
135 //! Standard reverse iterator.
136 typedef std::reverse_iterator<const_iterator> const_reverse_iterator;
138 /*! \brief
139 * Constructs an empty reference.
141 ArrayRef() : begin_(NULL), end_(NULL) {}
142 /*! \brief
143 * Constructs an empty reference.
145 * This is provided for convenience, such that EmptyArrayRef can be
146 * used to initialize any ArrayRef, without specifying the template
147 * type. It is not explicit to enable that usage.
149 ArrayRef(const EmptyArrayRef &) : begin_(nullptr), end_(nullptr) {}
150 /*! \brief
151 * Constructs a reference to a container or reference
153 * \param[in] o container to reference.
155 * Can be used to create a reference to a whole vector, std::array or
156 * an ArrayRef. The destination has to have a convertible pointer type
157 * (identical besides const or base class).
159 * Passed container must remain valid and not be reallocated for the
160 * lifetime of this object.
162 * This constructor is not explicit to allow directly passing
163 * a container to a method that takes ArrayRef.
165 template<typename U,
166 typename = typename std::enable_if<
167 std::is_convertible<typename std::remove_reference<U>::type::pointer,
168 pointer>::value>::type>
169 ArrayRef(U &&o) : begin_(o.data()), end_(o.data()+o.size()) {}
170 /*! \brief
171 * Constructs a reference to a particular range.
173 * \param[in] begin Pointer to the beginning of a range.
174 * \param[in] end Pointer to the end of a range.
176 * Passed pointers must remain valid for the lifetime of this object.
178 ArrayRef(pointer begin, pointer end)
179 : begin_(begin), end_(end)
181 GMX_ASSERT(end >= begin, "Invalid range");
183 //! \cond
184 // Doxygen 1.8.5 doesn't parse the declaration correctly...
185 /*! \brief
186 * Constructs a reference to a C array.
188 * \param[in] array C array to reference.
189 * \tparam count Deduced number of elements in \p array.
191 * This constructor can only be used with a real array (not with a
192 * pointer). It constructs a reference to the whole array, without
193 * a need to pass the number of elements explicitly. The compiler
194 * must be able to deduce the array size.
196 * Passed array must remain valid for the lifetime of this object.
198 * This constructor is not explicit to allow directly passing
199 * a C array to a function that takes an ArrayRef parameter.
201 template <size_t count>
202 ArrayRef(value_type (&array)[count])
203 : begin_(array), end_(array + count)
206 //! \endcond
208 //! Returns a reference to part of the container.
209 ArrayRef subArray(size_type start, size_type count) const
211 return {begin_+start, begin_+start+count};
213 //! Returns an iterator to the beginning of the container.
214 iterator begin() const { return begin_; }
215 //! Returns an iterator to the end of the container.
216 iterator end() const { return end_; }
217 //! Returns an iterator to the reverse beginning of the container.
218 reverse_iterator rbegin() const { return reverse_iterator(end()); }
219 //! Returns an iterator to the reverse end of the container.
220 reverse_iterator rend() const { return reverse_iterator(begin()); }
222 //! Returns the size of the container.
223 size_type size() const { return end_ - begin_; }
224 //! Identical to size().
225 size_type capacity() const { return end_ - begin_; }
226 //! Whether the container is empty.
227 bool empty() const { return begin_ == end_; }
229 //! Access container element.
230 reference operator[](size_type n) const { return begin_[n]; }
231 //! Access container element (throws on out-of-range error).
232 reference at(size_type n) const
234 if (n >= size())
236 throw std::out_of_range("Vector index out of range");
238 return begin_[n];
240 //! Returns the first element in the container.
241 reference front() const { return *begin_; }
242 //! Returns the first element in the container.
243 reference back() const { return *(end_ - 1); }
245 //! Returns a raw pointer to the contents of the array.
246 pointer data() const { return begin_; }
248 /*! \brief
249 * Swaps referenced memory with the other object.
251 * The actual memory areas are not modified, only the references are
252 * swapped.
254 void swap(ArrayRef<T> &other)
256 std::swap(begin_, other.begin_);
257 std::swap(end_, other.end_);
260 private:
261 pointer begin_;
262 pointer end_;
265 //! \copydoc ArrayRef::fromArray()
266 //! \related ArrayRef
267 template <typename T>
268 ArrayRef<T> arrayRefFromArray(T *begin, size_t size)
270 return ArrayRef<T>(begin, begin+size);
273 //! \copydoc ArrayRef::fromArray()
274 //! \related ArrayRef
275 template <typename T>
276 ArrayRef<const T> constArrayRefFromArray(const T *begin, size_t size)
278 return ArrayRef<const T>(begin, begin+size);
281 /*! \brief
282 * Create ArrayRef from container with type deduction
284 * \see ArrayRef
286 template <typename T>
287 ArrayRef<typename std::conditional<std::is_const<T>::value,
288 const typename T::value_type,
289 typename T::value_type>::type>
290 makeArrayRef(T &c)
292 return c;
295 /*! \brief
296 * Create ArrayRef to const T from container with type deduction
298 * \see ArrayRef
300 template <typename T>
301 ArrayRef<const typename T::value_type> makeConstArrayRef(const T &c)
303 return c;
306 /*! \brief
307 * Simple swap method for ArrayRef objects.
309 * \see ArrayRef::swap()
311 * \ingroup module_utility
313 template <typename T>
314 void swap(ArrayRef<T> &a, ArrayRef<T> &b)
316 a.swap(b);
319 /*! \brief Return a vector that is a copy of an ArrayRef.
321 * This makes it convenient, clear, and performant (the compiler will
322 * either do RVO to elide the temporary, or invoke the move constructor
323 * taking the unnamed temporary) to write a declaration like
325 * auto v = copyOf(arrayRef);
327 * \ingroup module_utility
329 template <typename T>
330 std::vector<T> copyOf(const ArrayRef<const T> &arrayRef)
332 return std::vector<T>(arrayRef.begin(), arrayRef.end());
335 } // namespace gmx
337 #endif