Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / math / multidimarray.h
blobe1e83fa924a6915f06550ecd154dcaa351ef6c7d
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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 /*! \libinternal
36 * \file
37 * \brief Declares MultiDimArray.
39 * \author Christian Blau <cblau@gwdg.de>
40 * \ingroup module_math
41 * \ingroup module_mdspan
44 #ifndef GMX_MATH_MULTIDIMARRAY_H_
45 #define GMX_MATH_MULTIDIMARRAY_H_
47 #include "gromacs/mdspan/mdspan.h"
48 #include "gromacs/utility/arrayref.h"
50 namespace gmx
53 namespace detail
55 //! Same as std::void_t from C++17
56 template< class ... >
57 using void_t = void;
59 template<typename T, typename = void>
60 struct is_resizable : std::false_type {};
62 template<typename T>
63 struct is_resizable<T, void_t<decltype(std::declval<T>().resize(size_t()))> > : std::true_type {};
65 //! Type has a resize member function callable with size_t argument
66 template<typename T>
67 constexpr bool is_resizable_v = is_resizable<T>::value;
68 } //namespace detail
70 /*! \libinternal \brief
71 * Multidimensional array that manages its own memory.
73 * \note No bounds checking when accessing memory
75 * \note That the view holds a valid pointer to the data is a class invariant.
77 * The Container type that stores the data may be resizable (std::vector or similar)
78 * or static (std::array or similar). Copy and move assignment routines as well as
79 * swapping are designed to yield good performances in both cases, notably
80 * foregoing the copy-swap idiom due to the bad performance in swapping std::array.
82 * This class avoids throwing exeptions, apart from the ones that might be thrown
83 * from the containers during resizing an allocation. (bad_malloc from std::vector
84 * is a likely candidate)
87 * Holds as many elements as required by a multidimensional view.
88 * \tparam TContainer Data type container for the data to be stored
89 * as MultiDimArray with random element access and
90 * value_type, refence and const_reference exposed
91 * \tparam Extents An extents class describing the array dimensions
92 * as used in module_mdspan
93 * \tparam LayoutPolicy The data layout as in module_mdspan describes
94 * translation of indices to memory offset.
95 * Defaults to right aligned, so that the right-most index
96 * is contiguous in memory.
98 template <class TContainer, class Extents, class LayoutPolicy = layout_right>
99 class MultiDimArray
101 public:
102 //! the type of values that are stored
103 using value_type = typename TContainer::value_type;
104 //! reference type to the stored values
105 using reference = typename TContainer::reference;
106 //! const reference type to the stored values
107 using const_reference = typename TContainer::const_reference;
108 //! the view used to access the data
109 using view_type = basic_mdspan<value_type, Extents, LayoutPolicy>;
110 //! const view on the data
111 using const_view_type = basic_mdspan<const value_type, Extents, LayoutPolicy>;
112 /*! \brief Iterator type for contiguous iteration over the stored data.
113 * Used, e.g., in free begin and end functions
115 using iterator = typename ArrayRef<value_type>::iterator;
116 /*! \brief Const iterator type for contiguous iteration over the stored data.
117 * used, e.g., in free begin and end functions
119 using const_iterator = const typename ArrayRef<const value_type>::const_iterator;
121 static_assert(detail::is_resizable_v<TContainer> == (Extents::rank_dynamic() > 0),
122 "Resizable container (e.g. std::vector) requires at least one dynamic rank. "
123 "Non-resizable container (e.g. std::array) requires zero dynamic ranks.");
125 /*! \brief
126 * Allocate dynamic array data and set view with the dynamic extents.
128 * \param[in] dynamicExtent A parameter pack that describes the dynamic
129 * size of the array. Empty if purely static.
131 * \tparam IndexType Parameter pack type holding the dynamic
132 * extents of the multidimensional array
134 template <class ... IndexType, typename T = TContainer,
135 typename = typename std::enable_if_t<detail::is_resizable_v<T> > >
136 MultiDimArray(IndexType... dynamicExtent)
138 resize(dynamicExtent ...);
140 /*! \brief
141 * Construction from fixed sized arrays if the array size is static and
142 * layout policy allows compile time determination of the container size.
144 * Enables the expected initialization
145 * MultiDimArray<std::array<float, 9>, extents<3,3>> arr = {{1,2...}}
146 * \tparam T Template parameter for activation via SFINAE.
148 // SFINAE required because std::vector::size isn't constexpr and is_constexpr doesn't exist.
149 template <typename T = TContainer,
150 typename = typename std::enable_if_t<!detail::is_resizable_v<T> > >
151 constexpr MultiDimArray(const TContainer &data = {}
152 ) noexcept : data_(data), view_(data_.data())
154 static_assert(TContainer().size() == typename view_type::mapping_type().required_span_size(),
155 "Non-resizable container type size must match static MultiDimArray size.");
157 //! Copy constructor
158 constexpr MultiDimArray(const MultiDimArray &o) :
159 data_(o.data_),
160 view_(data_.data(), o.view_.extents())
162 //! Move constructor
163 MultiDimArray(MultiDimArray &&o) noexcept : data_(std::move(o.data_)),
164 view_(data_.data(), o.view_.extents())
166 //! Copy assignment
167 MultiDimArray &operator=(const MultiDimArray &o)
169 data_ = o.data_;
170 view_ = view_type(data_.data(), o.view_.extents());
171 return *this;
173 //! Move assignment
174 MultiDimArray &operator=(MultiDimArray &&o) noexcept
176 data_ = std::move(o.data_);
177 view_ = view_type(data_.data(), o.view_.extents());
178 return *this;
180 //! Swaps content with other
181 void swap(MultiDimArray &o)
183 using std::swap;
184 swap(data_, o.data_);
185 // swap(view_, o.view_) also swaps the pointer to the data and thus does not work
186 // instead, restore the view as class invariant after the data swapping operation
187 o.view_ = view_type(o.data_.data(), view_.extents());
188 view_ = view_type(data_.data(), o.view_.extents());
190 /*! \brief
191 * Resize the dynamic extents of the array if any and set container size
192 * accordingly.
194 * Invalidates data and views of this array.
196 * \param[in] dynamicExtent A parameter pack that describes the dynamic
197 * size of the array. Empty if purely static.
198 * \tparam IndexType Parameter pack type holding the dynamic
199 * extents of the multidimensional array
201 template <class ... IndexType>
202 void resize(IndexType... dynamicExtent)
204 // use a mapping object to determine the required span size;
205 layout_right::mapping<Extents> map {Extents {dynamicExtent ...}};
206 data_.resize(map.required_span_size());
207 // to construct a valid view on the data, the container has to be resized before
208 // the assignment, so that data_.data() is valid
209 view_ = view_type(data_.data(), dynamicExtent ...);
211 /*! \brief Data access via multidimensional indices.
212 * This allows referencing rank R array elements as array(x_0,x_1,x_2, .., x_R)
214 * \param[in] index multidimensional indices as parameter pack
215 * the number of parameters must match the rank of the array.
217 * \returns reference to array element
219 template< class... IndexType>
220 reference operator()(IndexType ... index) noexcept
222 return view_(index ...);
224 /*! \brief Const data access via multidimensional indices.
225 * This allows referencing rank R array elements as array(x_0,x_1,x_2, .., x_R)
227 * \param[in] index multidimensional indices as parameter pack
228 * the number of parameters must match the rank of the array.
230 * \returns const reference to array element
232 template <class... IndexType>
233 constexpr const_reference operator()(IndexType... index) const noexcept
235 return view_(index ...);
237 /*! \brief Contiguous access to the data.
238 * \returns ArrayRef to stored data.
240 ArrayRef<value_type> toArrayRef()
242 return {data_.data(), data_.data()+data_.size()};
244 /*! \brief Contiguous const access to the data.
245 * \returns ArrayRef to stored data.
247 constexpr ArrayRef<const value_type> toArrayRef() const
249 return {data_.data(), data_.data()+data_.size()};
251 /*! \brief Return the extent.
252 * \param[in] k dimension to query for extent
253 * \returns extent along specified dimension
255 constexpr typename view_type::index_type extent(int k) const noexcept
257 return view_.extent(k);
259 //! Conversion to multidimensional view on the data
260 constexpr view_type asView() noexcept
262 return view_;
264 //! Conversion to const multidimensional view on the data
265 constexpr const_view_type asConstView() const noexcept
267 return {data_.data(), view_.mapping()};
269 private:
270 //! The contiguous data that is equipped with multidimensional indexing in this class
271 TContainer data_;
272 //! Multidimensional view into data_.
273 view_type view_;
276 //! Free MultiDimArray begin function addressing its contiguous memory.
277 template <class TContainer, class Extents>
278 constexpr typename MultiDimArray<TContainer, Extents>::const_iterator
279 begin(const MultiDimArray<TContainer, Extents> &multiDimArray)
281 return multiDimArray.toArrayRef().begin();
284 //! Free MultiDimArray begin function addressing its contiguous memory.
285 template <class TContainer, class Extents>
286 constexpr typename MultiDimArray<TContainer, Extents>::iterator
287 begin(MultiDimArray<TContainer, Extents> &multiDimArray)
289 return multiDimArray.toArrayRef().begin();
292 //! Free MultiDimArray end function addressing its contiguous memory.
293 template <class TContainer, class Extents>
294 constexpr typename MultiDimArray<TContainer, Extents>::const_iterator
295 end(const MultiDimArray<TContainer, Extents> &multiDimArray)
297 return multiDimArray.toArrayRef().end();
300 //! Free MultiDimArray end function addressing its contiguous memory.
301 template <class TContainer, class Extents>
302 constexpr typename MultiDimArray<TContainer, Extents>::iterator
303 end(MultiDimArray<TContainer, Extents> &multiDimArray)
305 return multiDimArray.toArrayRef().end();
308 //! Swap function
309 template <class TContainer, class Extents>
310 void swap(MultiDimArray<TContainer, Extents> &a, MultiDimArray<TContainer, Extents > &b) noexcept
312 a.swap(b);
315 } // namespace gmx
317 #endif // GMX_MATH_MULTIDIMARRAY_H_