Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / api / legacy / include / gromacs / utility / listoflists.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \file
36  * \brief
37  * Declares gmx::ListOfLists
38  *
39  * \author Berk Hess <hess@kth.se>
40  * \inpublicapi
41  * \ingroup module_utility
42  */
43 #ifndef GMX_UTILITY_LISTOFLISTS_H
44 #define GMX_UTILITY_LISTOFLISTS_H
45
46 #include <vector>
47
48 #include "gromacs/utility/arrayref.h"
49 #include "gromacs/utility/basedefinitions.h"
50 #include "gromacs/utility/exceptions.h"
51
52 namespace gmx
53 {
54
55 /*! \brief A list of lists, optimized for performance
56  *
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.
60  *
61  * This implementation stores all data internally in two std::vector objects
62  * and thereby avoids the overhead of managing \p size() separate objects
63  * in memory.
64  *
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_.
69  *
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.
73  *
74  * \tparam T value type
75  */
76
77 template<typename T>
78 class ListOfLists
79 {
80     // TODO: Use std::is_arithmetic_v when CUDA 11 is a requirement.
81     static_assert(std::is_arithmetic<T>::value, "This class is limited to arithmetic types");
82
83 public:
84     //! Constructs an empty list of lists
85     ListOfLists() = default;
86
87     /*! \brief Constructs a list of list from raw data in internal layout
88      *
89      * Does basic consistency checks and throws when one of those fail.
90      *
91      * \param[in] listRanges  Ranges of the lists concatenated (see above), is consumed
92      * \param[in] elements    Elements for all lists concatenated, is consumed
93      */
94     ListOfLists(std::vector<int>&& listRanges, std::vector<T>&& elements) :
95         listRanges_(std::move(listRanges)),
96         elements_(std::move(elements))
97     {
98         if (listRanges_.empty() || listRanges_.at(0) != 0)
99         {
100             GMX_THROW(InconsistentInputError(
101                     "listRanges does not have a first element with value 0"));
102         }
103         if (int(elements_.size()) != listRanges_.back())
104         {
105             GMX_THROW(InconsistentInputError(
106                     "The size of elements does not match the last value in listRanges"));
107         }
108     }
109
110     //! Returns the number of lists
111     std::size_t size() const { return listRanges_.size() - 1; }
112
113     /*! \brief Returns the number of lists
114      *
115      * \note Use ssize for any expression involving arithmetic operations
116      * (including loop indices).
117      */
118     index ssize() const { return index(listRanges_.size()) - 1; }
119
120     //! Returns whether the list holds no lists
121     bool empty() const { return listRanges_.size() == 1; }
122
123     //! Returns the sum of the number of elements over all lists
124     int numElements() const { return listRanges_.back(); }
125
126     //! Appends a new list with elements \p values, pass {} to add an empty list
127     void pushBack(ArrayRef<const T> values)
128     {
129         elements_.insert(elements_.end(), values.begin(), values.end());
130         listRanges_.push_back(int(elements_.size()));
131     }
132
133     //! Appends a new list with \p numElements elements
134     void pushBackListOfSize(int numElements)
135     {
136         // With arithmetic types enforced, this assertion is always true
137         // TODO: Use std::is_default_constructible_v when CUDA 11 is a requirement.
138         static_assert(std::is_default_constructible<T>::value,
139                       "pushBackListOfSize should only be called with default constructable types");
140         elements_.resize(elements_.size() + numElements);
141         listRanges_.push_back(int(elements_.size()));
142     }
143
144     //! Returns an ArrayRef to the elements of the list with the given index
145     ArrayRef<const T> operator[](std::size_t listIndex) const
146     {
147         return ArrayRef<const T>(elements_.data() + listRanges_[listIndex],
148                                  elements_.data() + listRanges_[listIndex + 1]);
149     }
150
151     //! Returns the list of elements for the list with index \p listIndex, throws an \p out_of_range exception when out of range
152     ArrayRef<const T> at(std::size_t listIndex) const
153     {
154         return ArrayRef<const T>(elements_.data() + listRanges_.at(listIndex),
155                                  elements_.data() + listRanges_.at(listIndex + 1));
156     }
157
158     /*! \brief Returns a reference to the first list
159      *
160      * \returns a reference to the first list
161      */
162     ArrayRef<T> front()
163     {
164         GMX_ASSERT(size() > 0, "Must contain a list if front() is called");
165         auto beginPtr = elements_.data();
166         auto endPtr   = beginPtr + listRanges_[1];
167         return { beginPtr, endPtr };
168     }
169     /*! \brief Returns a reference to the final list
170      *
171      * \returns a reference to the final list
172      */
173     ArrayRef<T> back()
174     {
175         GMX_ASSERT(size() > 0, "Must contain a list if bank() is called");
176         auto endIndex   = *(listRanges_.end() - 1);
177         auto beginIndex = *(listRanges_.end() - 2);
178         return { elements_.data() + beginIndex, elements_.data() + endIndex };
179     }
180
181     //! Clears the list
182     void clear()
183     {
184         listRanges_.resize(1);
185         elements_.clear();
186     }
187
188     //! Appends a ListOfLists at the end and increments the appended elements by \p offset
189     void appendListOfLists(const ListOfLists& listOfLists, const T offset = 0)
190     {
191         listRanges_.insert(
192                 listRanges_.end(), listOfLists.listRanges_.begin() + 1, listOfLists.listRanges_.end());
193         const int oldNumElements = elements_.size();
194         for (std::size_t i = listRanges_.size() - listOfLists.size(); i < listRanges_.size(); i++)
195         {
196             listRanges_[i] += oldNumElements;
197         }
198         elements_.insert(elements_.end(), listOfLists.elements_.begin(), listOfLists.elements_.end());
199
200         if (offset != 0)
201         {
202             for (std::size_t i = elements_.size() - listOfLists.elements_.size(); i < elements_.size(); i++)
203             {
204                 elements_[i] += offset;
205             }
206         }
207     }
208
209     //! Returns concatenated ranges of the lists (see above for details)
210     ArrayRef<const int> listRangesView() const { return listRanges_; }
211
212     //! Returns the a view of the elements of all lists concatenated
213     ArrayRef<const T> elementsView() const { return elements_; }
214
215 private:
216     //! The ranges of the lists, list i uses range \p listRanges_[i], \p listRanges_[i+1].
217     std::vector<int> listRanges_ = { 0 };
218     //! The elements in all lists concatenated
219     std::vector<T> elements_;
220 };
221
222 } // namespace gmx
223
224 #endif