--- /dev/null
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2019, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \file
+ * \brief
+ * Declares gmx::ListOfLists
+ *
+ * \author Berk Hess <hess@kth.se>
+ * \inpublicapi
+ * \ingroup module_utility
+ */
+#ifndef GMX_UTILITY_LISTOFLISTS_H
+#define GMX_UTILITY_LISTOFLISTS_H
+
+#include <vector>
+
+#include "gromacs/utility/arrayref.h"
+#include "gromacs/utility/basedefinitions.h"
+#include "gromacs/utility/exceptions.h"
+
+namespace gmx
+{
+
+/*! \brief A list of lists, optimized for performance
+ *
+ * This class holds a list of \p size() lists of elements of type \p T.
+ * To optimize performance, the only modification operation supporting
+ * is adding a new list at the end of the list of lists.
+ *
+ * This implementation stores all data internally in two std::vector objects
+ * and thereby avoids the overhead of managing \p size() separate objects
+ * in memory.
+ *
+ * Internal storage consists of one std::vector<int> listRanges_ of size number
+ * of lists plus one and a std::vector<T> elements_ with the elements of all
+ * lists concatenated. List i is stored in entries listRanges_[i] to
+ * listRanges_[i+1] in elements_.
+ *
+ * \tparam T value type
+ */
+template<typename T>
+class ListOfLists
+{
+public:
+ //! Constructs an empty list of lists
+ ListOfLists() = default;
+
+ /*! \brief Constructs a list of list from raw data in internal layout
+ *
+ * Does basic consistency checks and throws when one of those fail.
+ *
+ * \param[in] listRanges Ranges of the lists concatenated (see above), is consumed
+ * \param[in] elements Elements for all lists concatenated, is consumed
+ */
+ ListOfLists(std::vector<int>&& listRanges, std::vector<T>&& elements) :
+ listRanges_(std::move(listRanges)),
+ elements_(std::move(elements))
+ {
+ if (listRanges_.empty() || listRanges_.at(0) != 0)
+ {
+ GMX_THROW(InconsistentInputError(
+ "listRanges does not have a first element with value 0"));
+ }
+ if (int(elements_.size()) != listRanges_.back())
+ {
+ GMX_THROW(InconsistentInputError(
+ "The size of elements does not match the last value in listRanges"));
+ }
+ }
+
+ //! Returns the number of lists
+ std::size_t size() const { return listRanges_.size() - 1; }
+
+ /*! \brief Returns the number of lists
+ *
+ * \note Use ssize for any expression involving arithmetic operations
+ * (including loop indices).
+ */
+ index ssize() const { return index(listRanges_.size()) - 1; }
+
+ //! Returns whether the list holds no lists
+ bool empty() const { return listRanges_.size() == 1; }
+
+ //! Returns the sum of the number of elements over all lists
+ int numElements() const { return listRanges_.back(); }
+
+ //! Appends a new list with elements \p values, pass {} to add an empty list
+ void pushBack(ArrayRef<const T> values)
+ {
+ elements_.insert(elements_.end(), values.begin(), values.end());
+ listRanges_.push_back(int(elements_.size()));
+ }
+
+ //! Returns an ArrayRef to the elements of the list with the given index
+ ArrayRef<const T> operator[](std::size_t listIndex) const
+ {
+ return ArrayRef<const T>(elements_.data() + listRanges_[listIndex],
+ elements_.data() + listRanges_[listIndex + 1]);
+ }
+
+ //! Returns the list of elements for the list with index \p listIndex, throws an \p out_of_range exception when out of range
+ ArrayRef<const T> at(std::size_t listIndex) const
+ {
+ return ArrayRef<const T>(elements_.data() + listRanges_.at(listIndex),
+ elements_.data() + listRanges_.at(listIndex + 1));
+ }
+
+ //! Clears the list
+ void clear()
+ {
+ listRanges_.resize(1);
+ elements_.clear();
+ }
+
+ //! Appends a ListOfLists at the end
+ void appendListOfLists(const ListOfLists& listOfLists)
+ {
+ const std::size_t oldNumLists = size();
+ listRanges_.insert(listRanges_.end(), listOfLists.listRanges_.begin() + 1,
+ listOfLists.listRanges_.end());
+ const int oldNumElements = elements_.size();
+ for (std::size_t i = oldNumLists + 1; i < listRanges_.size(); i++)
+ {
+ listRanges_[i] += oldNumElements;
+ }
+ elements_.insert(elements_.end(), listOfLists.elements_.begin(), listOfLists.elements_.end());
+ }
+
+ //! Appends a ListOfLists at the end and increments the appended elements by \p offset
+ std::enable_if_t<std::is_arithmetic<T>::value, void> appendListOfLists(const ListOfLists& listOfLists,
+ const T offset)
+ {
+ const std::size_t oldNumElements = elements_.size();
+ appendListOfLists(listOfLists);
+ for (std::size_t i = oldNumElements; i < elements_.size(); i++)
+ {
+ elements_[i] += offset;
+ }
+ }
+
+ //! Returns concatenated ranges of the lists (see above for details)
+ ArrayRef<const int> listRangesView() const { return listRanges_; }
+
+ //! Returns the a view of the elements of all lists concatenated
+ ArrayRef<const T> elementsView() const { return elements_; }
+
+private:
+ //! The ranges of the lists, list i uses range \p listRanges_[i], \p listRanges_[i+1].
+ std::vector<int> listRanges_ = { 0 };
+ //! The elements in all lists concatenated
+ std::vector<T> elements_;
+};
+
+} // namespace gmx
+
+#endif
--- /dev/null
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2018,2019, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \internal \file
+ * \brief
+ * Tests for the ListOfLists class.
+ *
+ * \author berk Hess <hess@kth.se>
+ * \ingroup module_utility
+ */
+#include "gmxpre.h"
+
+#include "gromacs/utility/listoflists.h"
+
+#include <gtest/gtest.h>
+
+#include "testutils/testasserts.h"
+
+namespace gmx
+{
+
+namespace
+{
+
+//! Compares two lists element by element
+void compareLists(ArrayRef<const char> a, const std::vector<char>& b)
+{
+ EXPECT_EQ(a.size(), b.size());
+ if (a.size() == b.size())
+ {
+ for (size_t i = 0; i < b.size(); i++)
+ {
+ EXPECT_EQ(a[i], b[i]);
+ }
+ }
+}
+
+TEST(ListOfLists, EmptyListOfListsWorks)
+{
+ ListOfLists<char> list;
+
+ EXPECT_EQ(list.size(), 0);
+ EXPECT_EQ(list.empty(), true);
+ EXPECT_EQ(list.numElements(), 0);
+}
+
+TEST(ListOfLists, AppendWorks)
+{
+ ListOfLists<char> list;
+
+ std::vector<char> v1 = { 5, 3 };
+ std::vector<char> v2 = { -1, 7, 4 };
+ list.pushBack(v1);
+ list.pushBack(v2);
+ EXPECT_EQ(list.size(), 2);
+ compareLists(list[0], v1);
+ compareLists(list[1], v2);
+}
+
+TEST(ListOfLists, EmptyListWorks)
+{
+ ListOfLists<char> list;
+
+ std::vector<char> v = { 5, 3 };
+ list.pushBack(v);
+ list.pushBack({});
+ EXPECT_EQ(list.size(), 2);
+ auto a = list[1];
+ EXPECT_EQ(a.empty(), true);
+}
+
+TEST(ListOfLists, ClearWorks)
+{
+ ListOfLists<char> list;
+
+ std::vector<char> v = { 5, 3 };
+ list.pushBack(v);
+ list.pushBack({});
+ list.clear();
+ EXPECT_EQ(list.empty(), true);
+ EXPECT_EQ(list.numElements(), 0);
+}
+
+TEST(ListOfLists, OutOfRangeAccessThrows)
+{
+ ListOfLists<char> list;
+
+ std::vector<char> v = { 5, 3 };
+ EXPECT_THROW(list.at(1), std::out_of_range);
+}
+
+TEST(ListOfLists, ExtractsAndRestores)
+{
+ ListOfLists<char> list1;
+
+ std::vector<char> v1 = { 5, 3 };
+ std::vector<char> v3 = { -1, 4 };
+ list1.pushBack(v1);
+ list1.pushBack({});
+ list1.pushBack(v3);
+ auto listRanges = list1.listRangesView();
+ auto elements = list1.elementsView();
+ std::vector<int> listRangesVector;
+ listRangesVector.insert(listRangesVector.begin(), listRanges.begin(), listRanges.end());
+ std::vector<char> elementsVector;
+ elementsVector.insert(elementsVector.begin(), elements.begin(), elements.end());
+ ListOfLists<char> list2(std::move(listRangesVector), std::move(elementsVector));
+ compareLists(list2[0], v1);
+ EXPECT_EQ(list2[1].empty(), true);
+ compareLists(list2[2], v3);
+}
+
+// Test that we can extract raw data from one object and use it correctly generate a new object
+TEST(ListOfLists, AppendsListOfLists)
+{
+ ListOfLists<char> list1;
+ ListOfLists<char> list2;
+
+ std::vector<char> v1 = { 5, 3 };
+ list1.pushBack(v1);
+ std::vector<char> v2 = { 2, -1 };
+ std::vector<char> v3 = { 4 };
+ list2.pushBack(v2);
+ list2.pushBack(v3);
+ const char offset = 2;
+ list1.appendListOfLists(list2, offset);
+ EXPECT_EQ(list1.size(), 3);
+ auto a = list1[1];
+ EXPECT_EQ(a.size(), 2);
+ EXPECT_EQ(a[0], v2[0] + offset);
+ EXPECT_EQ(a[1], v2[1] + offset);
+}
+
+} // namespace
+
+} // namespace gmx