Specialize ArrayRef for SimdReal
[gromacs.git] / src / gromacs / utility / arrayref.h
blobb82c057ac482a30ad4f98de6dec15a599f212e9c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017, 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 private:
115 typedef typename std::remove_const<T>::type non_const_value_type;
116 public:
117 //! Type of values stored in the container.
118 typedef T value_type;
119 //! Type for representing size of the container.
120 typedef size_t size_type;
121 //! Type for representing difference between two container indices.
122 typedef ptrdiff_t difference_type;
123 //! Const reference to a container element.
124 typedef const T &const_reference;
125 //! Const pointer to a container element.
126 typedef const T *const_pointer;
127 //! Const iterator type for the container.
128 typedef const T *const_iterator;
129 //! Reference to a container element.
130 typedef T &reference;
131 //! Pointer to a container element.
132 typedef T *pointer;
133 //! Iterator type for the container.
134 typedef T *iterator;
135 //! Standard reverse iterator.
136 typedef std::reverse_iterator<iterator> reverse_iterator;
137 //! Standard reverse iterator.
138 typedef std::reverse_iterator<const_iterator> const_reverse_iterator;
140 /*! \brief
141 * Constructs a reference to a particular range from two pointers.
143 * \param[in] begin Pointer to the beginning of a range.
144 * \param[in] end Pointer to the end of a range.
146 * Passed pointers must remain valid for the lifetime of this object.
148 static ArrayRef<value_type>
149 fromPointers(value_type *begin, value_type *end)
151 return ArrayRef<value_type>(begin, end);
153 /*! \brief
154 * Constructs a reference to an array.
156 * \param[in] begin Pointer to the beginning of the array.
157 * May be NULL if \p size is zero.
158 * \param[in] size Number of elements in the array.
160 * Passed pointer must remain valid for the lifetime of this object.
162 static ArrayRef<value_type>
163 fromArray(value_type *begin, size_t size)
165 return ArrayRef<value_type>(begin, begin+size);
167 /*! \brief
168 * Constructs a reference to a particular range in a std::vector.
170 * \param[in] begin Iterator to the beginning of a range.
171 * \param[in] end Iterator to the end of a range.
173 * The referenced vector must remain valid and not be reallocated for
174 * the lifetime of this object.
176 static ArrayRef<value_type>
177 fromVector(typename std::vector<non_const_value_type>::iterator begin,
178 typename std::vector<non_const_value_type>::iterator end)
180 value_type *p_begin = (begin != end) ? &*begin : nullptr;
181 value_type *p_end = p_begin + (end-begin);
182 return ArrayRef<value_type>(p_begin, p_end);
184 //! \copydoc ArrayRef::fromVector(typename std::vector<non_const_value_type>::iterator, typename std::vector<non_const_value_type>::iterator)
185 static ArrayRef<value_type>
186 fromVector(typename std::vector<non_const_value_type>::const_iterator begin,
187 typename std::vector<non_const_value_type>::const_iterator end)
189 value_type *p_begin = (begin != end) ? &*begin : nullptr;
190 value_type *p_end = p_begin + (end-begin);
191 return ArrayRef<value_type>(p_begin, p_end);
193 /*! \brief
194 * Constructs an empty reference.
196 ArrayRef() : begin_(NULL), end_(NULL) {}
197 /*! \brief
198 * Constructs an empty reference.
200 * This is provided for convenience, such that EmptyArrayRef can be
201 * used to initialize any ArrayRef, without specifying the template
202 * type. It is not explicit to enable that usage.
204 ArrayRef(const EmptyArrayRef &) : begin_(nullptr), end_(nullptr) {}
205 /*! \brief
206 * Constructs a reference to const data from a reference to non-const data.
208 * Constructs a ArrayRef<const T> from a ArrayRef<T>.
210 template<typename = T> //Otherwise useless template argument
211 //to avoid this being used as copy constructor
212 ArrayRef(const ArrayRef<non_const_value_type> &o) :
213 begin_(o.begin()), end_(o.end()) {}
214 /*! \brief
215 * Constructs a reference to a particular range.
217 * \param[in] begin Pointer to the beginning of a range.
218 * \param[in] end Pointer to the end of a range.
220 * Passed pointers must remain valid for the lifetime of this object.
222 * \note For clarity, use the non-member function arrayRefFromPointers
223 * instead.
225 ArrayRef(pointer begin, pointer end)
226 : begin_(begin), end_(end)
228 GMX_ASSERT(end >= begin, "Invalid range");
230 /*! \brief
231 * Constructs a reference to a whole std::vector<T, A>.
233 * \param[in] v Vector to reference.
235 * Passed vector must remain valid and not be reallocated for the
236 * lifetime of this object.
238 * This constructor is not explicit to allow directly passing
239 * std::vector<T, A> to a method that takes ArrayRef.
241 template <typename A>
242 ArrayRef(std::vector<non_const_value_type, A> &v)
243 : begin_((!v.empty()) ? &v[0] : nullptr),
244 end_((!v.empty()) ? &v[0] + v.size() : nullptr)
247 //! \copydoc ArrayRef::ArrayRef(std::vector<non_const_value_type, A>&)
248 template <typename A>
249 ArrayRef(const std::vector<non_const_value_type, A> &v)
250 : begin_((!v.empty()) ? &v[0] : nullptr),
251 end_((!v.empty()) ? &v[0] + v.size() : nullptr)
254 /*! \brief
255 * Constructs a reference to a whole std::array<T>.
257 * \param[in] a Array to reference.
259 * Passed array must remain valid for the lifetime of this
260 * object.
262 * This constructor is not explicit to allow directly passing
263 * std::array<T> to a method that takes ArrayRef.
265 template <size_t count>
266 ArrayRef(std::array<non_const_value_type, count> &a)
267 : begin_((!a.empty()) ? &a[0] : NULL),
268 end_((!a.empty()) ? &a[0] + a.size() : NULL)
271 //! \copydoc ArrayRef::ArrayRef(std::array<non_const_value_type, count> &)
272 template <size_t count>
273 ArrayRef(const std::array<non_const_value_type, count> &a)
274 : begin_((!a.empty()) ? &a[0] : NULL),
275 end_((!a.empty()) ? &a[0] + a.size() : NULL)
278 //! \cond
279 // Doxygen 1.8.5 doesn't parse the declaration correctly...
280 /*! \brief
281 * Constructs a reference to a C array.
283 * \param[in] array C array to reference.
284 * \tparam count Deduced number of elements in \p array.
286 * This constructor can only be used with a real array (not with a
287 * pointer). It constructs a reference to the whole array, without
288 * a need to pass the number of elements explicitly. The compiler
289 * must be able to deduce the array size.
291 * Passed array must remain valid for the lifetime of this object.
293 * This constructor is not explicit to allow directly passing
294 * a C array to a function that takes an ArrayRef parameter.
296 template <size_t count>
297 ArrayRef(value_type (&array)[count])
298 : begin_(array), end_(array + count)
301 //! \endcond
303 //! Returns an iterator to the beginning of the container.
304 iterator begin() const { return begin_; }
305 //! Returns an iterator to the end of the container.
306 iterator end() const { return end_; }
307 //! Returns an iterator to the reverse beginning of the container.
308 reverse_iterator rbegin() const { return reverse_iterator(end()); }
309 //! Returns an iterator to the reverse end of the container.
310 reverse_iterator rend() const { return reverse_iterator(begin()); }
312 //! Returns the size of the container.
313 size_type size() const { return end_ - begin_; }
314 //! Identical to size().
315 size_type capacity() const { return end_ - begin_; }
316 //! Whether the container is empty.
317 bool empty() const { return begin_ == end_; }
319 //! Access container element.
320 reference operator[](size_type n) const { return begin_[n]; }
321 //! Access container element (throws on out-of-range error).
322 reference at(size_type n) const
324 if (n >= size())
326 throw std::out_of_range("Vector index out of range");
328 return begin_[n];
330 //! Returns the first element in the container.
331 reference front() const { return *begin_; }
332 //! Returns the first element in the container.
333 reference back() const { return *(end_ - 1); }
335 //! Returns a raw pointer to the contents of the array.
336 pointer data() const { return begin_; }
338 /*! \brief
339 * Swaps referenced memory with the other object.
341 * The actual memory areas are not modified, only the references are
342 * swapped.
344 void swap(ArrayRef<T> &other)
346 std::swap(begin_, other.begin_);
347 std::swap(end_, other.end_);
350 private:
351 pointer begin_;
352 pointer end_;
356 //! Alias for ArrayRef<const T>
357 //! \related ArrayRef
358 template <typename T>
359 using ConstArrayRef = ArrayRef<const T>;
362 //! \copydoc ArrayRef::fromPointers()
363 //! \related ArrayRef
364 template <typename T>
365 ArrayRef<T> arrayRefFromPointers(T *begin, T *end)
367 return ArrayRef<T>::fromPointers(begin, end);
369 //! \copydoc ArrayRef::fromArray()
370 //! \related ArrayRef
371 template <typename T>
372 ArrayRef<T> arrayRefFromArray(T *begin, size_t size)
374 return ArrayRef<T>::fromArray(begin, size);
376 //! \copydoc ArrayRef::fromVector()
377 //! \related ArrayRef
378 template <typename T>
379 ArrayRef<T> arrayRefFromVector(typename std::vector<T>::iterator begin,
380 typename std::vector<T>::iterator end)
382 return ArrayRef<T>::fromVector(begin, end);
386 //! \copydoc ArrayRef::fromPointers()
387 //! \related ArrayRef
388 template <typename T>
389 ConstArrayRef<T> constArrayRefFromPointers(const T *begin, const T *end)
391 return ConstArrayRef<T>::fromPointers(begin, end);
393 //! \copydoc ArrayRef::fromArray()
394 //! \related ArrayRef
395 template <typename T>
396 ConstArrayRef<T> constArrayRefFromArray(const T *begin, size_t size)
398 return ConstArrayRef<T>::fromArray(begin, size);
400 //! \copydoc ArrayRef::fromVector()
401 //! \related ArrayRef
402 template <typename T>
403 ConstArrayRef<T> constArrayRefFromVector(typename std::vector<T>::const_iterator begin,
404 typename std::vector<T>::const_iterator end)
406 return ConstArrayRef<T>::fromVector(begin, end);
409 /*! \brief
410 * Simple swap method for ArrayRef objects.
412 * \see ArrayRef::swap()
414 * \ingroup module_utility
416 template <typename T>
417 void swap(ArrayRef<T> &a, ArrayRef<T> &b)
419 a.swap(b);
422 } // namespace gmx
424 #endif