Modernize syntax of enable_if and traits to use _t helpers
[gromacs.git] / src / gromacs / utility / arrayref.h
blobac7a8a0c7dd3f21bc22253fa30563a5b157b34c9
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, 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 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 * \author Roland Schulz <roland.schulz@intel.com>
42 * \author Berk Hess <hess@kth.se>
43 * \inpublicapi
44 * \ingroup module_utility
46 #ifndef GMX_UTILITY_ARRAYREF_H
47 #define GMX_UTILITY_ARRAYREF_H
49 #include <cstddef>
51 #include <array>
52 #include <iterator>
53 #include <stdexcept>
54 #include <utility>
55 #include <vector>
57 #include "gromacs/utility/gmxassert.h"
59 namespace gmx
62 /*! \brief STL-like interface to a C array of T (or part
63 * of a std container of T).
65 * \tparam T Value type of elements.
67 * This class provides an interface similar to \c std::vector<T, A>, with the
68 * following main differences:
69 * - This class does not have its own storage. Instead, it references an
70 * existing array of values (either a C-style array or part of an existing
71 * std::vector<T, A> or std::array<T>).
72 * - It is only possible to modify the values themselves through ArrayRef;
73 * it is not possible to add or remove values.
74 * - Copying objects of this type is cheap, and the copies behave identically
75 * to the original object: the copy references the same set of values.
77 * This class is useful for writing wrappers that expose a view of the
78 * internal data stored as a single vector/array, which can be a whole
79 * or part of the underlying storage.
81 * Methods in this class do not throw, except where indicated.
83 * Note that due to a Doxygen limitation, the constructor that takes a C array
84 * whose size is known at compile time does not appear in the documentation.
86 * To refer to const data of type T, ArrayRef<const T> is used. For both const
87 * and non-const std::vector and std::array an ArrayRef view can be created.
88 * Attempting to create a non-const ArrayRef of a const vector/array will result
89 * in a compiler error in the respective constructor.
91 * For SIMD types there is template specialization available
92 * (e.g. ArrayRef<SimdReal>) in gromacs/simd/simd_memory.h which should have
93 * the same functionality as much as possible.
95 * \todo
96 * This class is not complete. There are likely also methods missing (not
97 * required for current usage).
99 * \inpublicapi
100 * \ingroup module_utility
102 template <typename T>
103 class ArrayRef
105 public:
106 //! Type of values stored in the reference.
107 typedef T value_type;
108 //! Type for representing size of the reference.
109 typedef size_t size_type;
110 //! Type for representing difference between two indices.
111 typedef ptrdiff_t difference_type;
112 //! Const reference to an element.
113 typedef const T &const_reference;
114 //! Const pointer to an element.
115 typedef const T *const_pointer;
116 //! Const iterator type to an element.
117 typedef const T *const_iterator;
118 //! Reference to an element.
119 typedef T &reference;
120 //! Pointer to an element.
121 typedef T *pointer;
122 //! Iterator type to an element.
123 typedef T *iterator;
124 //! Standard reverse iterator.
125 typedef std::reverse_iterator<iterator> reverse_iterator;
126 //! Standard reverse iterator.
127 typedef std::reverse_iterator<const_iterator> const_reverse_iterator;
129 /*! \brief
130 * Constructs an empty reference.
132 ArrayRef() : begin_(nullptr), end_(nullptr) {}
133 /*! \brief
134 * Constructs a reference to a container or reference
136 * \param[in] o container to reference.
138 * Can be used to create a reference to a whole vector, std::array or
139 * an ArrayRef. The destination has to have a convertible pointer type
140 * (identical besides const or base class).
142 * Passed container must remain valid and not be reallocated for the
143 * lifetime of this object.
145 * This constructor is not explicit to allow directly passing
146 * a container to a method that takes ArrayRef.
148 template<typename U,
149 typename = std::enable_if_t<
150 std::is_convertible<typename std::remove_reference_t<U>::pointer,
151 pointer>::value> >
152 ArrayRef(U &&o) : begin_(o.data()), end_(o.data()+o.size()) {}
153 /*! \brief
154 * Constructs a reference to a particular range.
156 * \param[in] begin Pointer to the beginning of a range.
157 * \param[in] end Pointer to the end of a range.
159 * Passed pointers must remain valid for the lifetime of this object.
161 ArrayRef(pointer begin, pointer end)
162 : begin_(begin), end_(end)
164 GMX_ASSERT(end >= begin, "Invalid range");
166 //! \cond
167 // Doxygen 1.8.5 doesn't parse the declaration correctly...
168 /*! \brief
169 * Constructs a reference to a C array.
171 * \param[in] array C array to reference.
172 * \tparam count Deduced number of elements in \p array.
174 * This constructor can only be used with a real array (not with a
175 * pointer). It constructs a reference to the whole array, without
176 * a need to pass the number of elements explicitly. The compiler
177 * must be able to deduce the array size.
179 * Passed array must remain valid for the lifetime of this object.
181 * This constructor is not explicit to allow directly passing
182 * a C array to a function that takes an ArrayRef parameter.
184 template <size_t count>
185 ArrayRef(value_type (&array)[count])
186 : begin_(array), end_(array + count)
189 //! \endcond
191 //! Returns a reference to part of the memory.
192 ArrayRef subArray(size_type start, size_type count) const
194 return {begin_+start, begin_+start+count};
196 //! Returns an iterator to the beginning of the reference.
197 iterator begin() const { return begin_; }
198 //! Returns an iterator to the end of the reference.
199 iterator end() const { return end_; }
200 //! Returns an iterator to the reverse beginning of the reference.
201 reverse_iterator rbegin() const { return reverse_iterator(end()); }
202 //! Returns an iterator to the reverse end of the reference.
203 reverse_iterator rend() const { return reverse_iterator(begin()); }
205 /*! \brief Returns the size of the reference.
207 * \note Use ssize for any expression involving arithmetic operations
208 (including loop indices).
210 size_type size() const { return end_ - begin_; }
211 //! Returns the signed size of the reference.
212 index ssize() const { return size(); }
213 //! Identical to size().
214 size_type capacity() const { return end_ - begin_; }
215 //! Whether the reference refers to no memory.
216 bool empty() const { return begin_ == end_; }
218 //! Access an element.
219 reference operator[](size_type n) const { return begin_[n]; }
220 //! Access an element (throws on out-of-range error).
221 reference at(size_type n) const
223 if (n >= size())
225 throw std::out_of_range("Vector index out of range");
227 return begin_[n];
229 //! Returns the first element.
230 reference front() const { return *begin_; }
231 //! Returns the first element.
232 reference back() const { return *(end_ - 1); }
234 //! Returns a raw pointer to the contents of the array.
235 pointer data() const { return begin_; }
237 /*! \brief
238 * Swaps referenced memory with the other object.
240 * The actual memory areas are not modified, only the references are
241 * swapped.
243 void swap(ArrayRef<T> &other)
245 std::swap(begin_, other.begin_);
246 std::swap(end_, other.end_);
249 private:
250 pointer begin_;
251 pointer end_;
254 //! \copydoc ArrayRef::fromArray()
255 //! \related ArrayRef
256 template <typename T>
257 ArrayRef<T> arrayRefFromArray(T *begin, size_t size)
259 return ArrayRef<T>(begin, begin+size);
262 //! \copydoc ArrayRef::fromArray()
263 //! \related ArrayRef
264 template <typename T>
265 ArrayRef<const T> constArrayRefFromArray(const T *begin, size_t size)
267 return ArrayRef<const T>(begin, begin+size);
270 /*! \brief
271 * Create ArrayRef from container with type deduction
273 * \see ArrayRef
275 template <typename T>
276 ArrayRef < std::conditional_t < std::is_const<T>::value,
277 const typename T::value_type,
278 typename T::value_type>>
279 makeArrayRef(T &c)
281 return c;
284 /*! \brief
285 * Create ArrayRef to const T from container with type deduction
287 * \see ArrayRef
289 template <typename T>
290 ArrayRef<const typename T::value_type> makeConstArrayRef(const T &c)
292 return c;
295 /*! \brief
296 * Simple swap method for ArrayRef objects.
298 * \see ArrayRef::swap()
300 * \ingroup module_utility
302 template <typename T>
303 void swap(ArrayRef<T> &a, ArrayRef<T> &b)
305 a.swap(b);
308 /*! \brief Return a vector that is a copy of an ArrayRef.
310 * This makes it convenient, clear, and performant (the compiler will
311 * either do RVO to elide the temporary, or invoke the move constructor
312 * taking the unnamed temporary) to write a declaration like
314 * auto v = copyOf(arrayRef);
316 * \ingroup module_utility
318 template <typename T>
319 std::vector<T> copyOf(const ArrayRef<const T> &arrayRef)
321 return std::vector<T>(arrayRef.begin(), arrayRef.end());
324 } // namespace gmx
326 #endif