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.
37 * Declares gmx::ListOfLists
39 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_utility
43 #ifndef GMX_UTILITY_LISTOFLISTS_H
44 #define GMX_UTILITY_LISTOFLISTS_H
48 #include "gromacs/utility/arrayref.h"
49 #include "gromacs/utility/basedefinitions.h"
50 #include "gromacs/utility/exceptions.h"
55 /*! \brief A list of lists, optimized for performance
57 * This class holds a list of \p size() lists of elements of type \p T.
58 * To optimize performance, the only modification operation supporting
59 * is adding a new list at the end of the list of lists.
61 * This implementation stores all data internally in two std::vector objects
62 * and thereby avoids the overhead of managing \p size() separate objects
65 * Internal storage consists of one std::vector<int> listRanges_ of size number
66 * of lists plus one and a std::vector<T> elements_ with the elements of all
67 * lists concatenated. List i is stored in entries listRanges_[i] to
68 * listRanges_[i+1] in elements_.
70 * \note This class is currently limited to arithmetic types, mainly because
71 * this should only be used for performance critical applications.
72 * When performance is not critical, a std::vector of std::vector can be used.
74 * \tparam T value type
80 static_assert(std::is_arithmetic
<T
>::value
, "This class is limited to arithmetic types");
83 //! Constructs an empty list of lists
84 ListOfLists() = default;
86 /*! \brief Constructs a list of list from raw data in internal layout
88 * Does basic consistency checks and throws when one of those fail.
90 * \param[in] listRanges Ranges of the lists concatenated (see above), is consumed
91 * \param[in] elements Elements for all lists concatenated, is consumed
93 ListOfLists(std::vector
<int>&& listRanges
, std::vector
<T
>&& elements
) :
94 listRanges_(std::move(listRanges
)),
95 elements_(std::move(elements
))
97 if (listRanges_
.empty() || listRanges_
.at(0) != 0)
99 GMX_THROW(InconsistentInputError(
100 "listRanges does not have a first element with value 0"));
102 if (int(elements_
.size()) != listRanges_
.back())
104 GMX_THROW(InconsistentInputError(
105 "The size of elements does not match the last value in listRanges"));
109 //! Returns the number of lists
110 std::size_t size() const { return listRanges_
.size() - 1; }
112 /*! \brief Returns the number of lists
114 * \note Use ssize for any expression involving arithmetic operations
115 * (including loop indices).
117 index
ssize() const { return index(listRanges_
.size()) - 1; }
119 //! Returns whether the list holds no lists
120 bool empty() const { return listRanges_
.size() == 1; }
122 //! Returns the sum of the number of elements over all lists
123 int numElements() const { return listRanges_
.back(); }
125 //! Appends a new list with elements \p values, pass {} to add an empty list
126 void pushBack(ArrayRef
<const T
> values
)
128 elements_
.insert(elements_
.end(), values
.begin(), values
.end());
129 listRanges_
.push_back(int(elements_
.size()));
132 //! Appends a new list with \p numElements elements
133 void pushBackListOfSize(int numElements
)
135 // With arithmetic types enforced, this assertion is always true
136 static_assert(std::is_default_constructible
<T
>::value
,
137 "pushBackListOfSize should only be called with default constructable types");
138 elements_
.resize(elements_
.size() + numElements
);
139 listRanges_
.push_back(int(elements_
.size()));
142 //! Returns an ArrayRef to the elements of the list with the given index
143 ArrayRef
<const T
> operator[](std::size_t listIndex
) const
145 return ArrayRef
<const T
>(elements_
.data() + listRanges_
[listIndex
],
146 elements_
.data() + listRanges_
[listIndex
+ 1]);
149 //! Returns the list of elements for the list with index \p listIndex, throws an \p out_of_range exception when out of range
150 ArrayRef
<const T
> at(std::size_t listIndex
) const
152 return ArrayRef
<const T
>(elements_
.data() + listRanges_
.at(listIndex
),
153 elements_
.data() + listRanges_
.at(listIndex
+ 1));
156 /*! \brief Returns a reference to the first list
158 * \returns a reference to the first list
162 GMX_ASSERT(size() > 0, "Must contain a list if front() is called");
163 auto beginPtr
= elements_
.data();
164 auto endPtr
= beginPtr
+ listRanges_
[1];
165 return { beginPtr
, endPtr
};
167 /*! \brief Returns a reference to the final list
169 * \returns a reference to the final list
173 GMX_ASSERT(size() > 0, "Must contain a list if bank() is called");
174 auto endIndex
= *(listRanges_
.end() - 1);
175 auto beginIndex
= *(listRanges_
.end() - 2);
176 return { elements_
.data() + beginIndex
, elements_
.data() + endIndex
};
182 listRanges_
.resize(1);
186 //! Appends a ListOfLists at the end and increments the appended elements by \p offset
187 void appendListOfLists(const ListOfLists
& listOfLists
, const T offset
= 0)
189 listRanges_
.insert(listRanges_
.end(), listOfLists
.listRanges_
.begin() + 1,
190 listOfLists
.listRanges_
.end());
191 const int oldNumElements
= elements_
.size();
192 for (std::size_t i
= listRanges_
.size() - listOfLists
.size(); i
< listRanges_
.size(); i
++)
194 listRanges_
[i
] += oldNumElements
;
196 elements_
.insert(elements_
.end(), listOfLists
.elements_
.begin(), listOfLists
.elements_
.end());
200 for (std::size_t i
= elements_
.size() - listOfLists
.elements_
.size(); i
< elements_
.size(); i
++)
202 elements_
[i
] += offset
;
207 //! Returns concatenated ranges of the lists (see above for details)
208 ArrayRef
<const int> listRangesView() const { return listRanges_
; }
210 //! Returns the a view of the elements of all lists concatenated
211 ArrayRef
<const T
> elementsView() const { return elements_
; }
214 //! The ranges of the lists, list i uses range \p listRanges_[i], \p listRanges_[i+1].
215 std::vector
<int> listRanges_
= { 0 };
216 //! The elements in all lists concatenated
217 std::vector
<T
> elements_
;