From: David Wells Date: Sun, 18 Mar 2018 22:41:09 +0000 (-0400) Subject: Implement a random-access iterator for TransposeTable. X-Git-Tag: v9.0.0-rc1~290^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=b4b71dc9f60;p=dealii.git Implement a random-access iterator for TransposeTable. --- diff --git a/doc/news/changes/minor/20180225DavidWells b/doc/news/changes/minor/20180225DavidWells new file mode 100644 index 0000000000..503a5fe89a --- /dev/null +++ b/doc/news/changes/minor/20180225DavidWells @@ -0,0 +1,4 @@ +New: TransposeTable (the base class of LAPACKFullMatrix and ScaLAPACKMatrix) now +has a random access iterator implementation similar to the one provided by +SparseMatrix. +
(David Wells, 2018/02/25) diff --git a/include/deal.II/base/table.h b/include/deal.II/base/table.h index f4c8b05bae..cf263b1059 100644 --- a/include/deal.II/base/table.h +++ b/include/deal.II/base/table.h @@ -22,6 +22,7 @@ #include #include #include +#include #include #include @@ -31,6 +32,7 @@ DEAL_II_NAMESPACE_OPEN // forward declaration template class TableBase; template class Table; +template class TransposeTable; template class Table<1,T>; template class Table<2,T>; template class Table<3,T>; @@ -1499,6 +1501,202 @@ public: +/** + * A namespace for iterators for TransposeTable. TransposeTable is unique in + * that it stores entries in column-major order. + */ +namespace TransposeTableIterators +{ + // Forward declaration of the iterator class. + template + class Iterator; + + /** + * Accessor class template. This will be partially specialized for both + * values of Constness. + */ + template + class Accessor; + + /** + * Accessor base class. Stores a const pointer when Constness == + * true and a non-const pointer when Constness == false. + */ + template + class AccessorBase + { + public: + /** + * Type of the stored pointer to the TransposeTable. + */ + typedef typename std::conditional*, TransposeTable*>::type + container_pointer_type; + + /** + * Default constructor. + */ + AccessorBase(); + + /** + * Constructor setting up the end iterator. + */ + AccessorBase(const container_pointer_type table); + + /** + * Copy constructor from a non-const Accessor. + */ + AccessorBase(const AccessorBase &); + + /** + * Constructor taking an array index. Note that TransposeTable stores data + * in column-major order. + */ + AccessorBase(const container_pointer_type table, + const std::ptrdiff_t linear_index); + + /** + * Get a constant reference to the value of the element represented by + * this object. + */ + const T &value() const; + + /** + * Numerical type of the row and column indices of the TransposeTable. + */ + typedef typename TransposeTable::size_type size_type; + + /** + * Access the row of the current entry. + */ + size_type row() const; + + /** + * Access the column of the current entry. + */ + size_type column() const; + + protected: + /** + * Pointer to the TransposeTable. + */ + container_pointer_type container; + + /** + * Current index. + */ + std::ptrdiff_t linear_index; + + /** + * Make the const version a friend for copying. + */ + friend class AccessorBase; + + /** + * Make the underlying iterator class a friend. + */ + friend class LinearIndexIterator, Accessor>; + }; + + /** + * Constant Accessor. Same as the base class. + */ + template + class Accessor : public AccessorBase + { + public: + /** + * Inherit the base class constructors. + */ + using AccessorBase::AccessorBase; + }; + + /** + * Accessor class for non-constant tables. Defines a few extra operations + * for modifying values. + */ + template + class Accessor : public AccessorBase + { + public: + /** + * Inherit the base class constructors. + */ + using AccessorBase::AccessorBase; + + /** + * Assignment operator. This assigns a new value to the table entry at the + * current row and column coordinates. + */ + const Accessor &operator = (const T &) const; + + /** + * Move assignment operator. This assigns a new value to the table entry at the + * current row and column coordinates. + */ + const Accessor &operator = (T&&) const; + + /** + * Since we overload value() we have to explicitly use the base class + * version. + */ + using AccessorBase::value; + + /** + * Get a reference to the value of the element represented by + * this object. + */ + T &value() const; + }; + + /** + * Iterator class for TransposeTable. + */ + template + class Iterator : public LinearIndexIterator, Accessor> + { + public: + /** + * Size type used by TransposeTable. + */ + typedef typename TransposeTable::size_type size_type; + + /** + * Type of the stored pointer to the TransposeTable. + */ + typedef typename std::conditional*, TransposeTable*>::type + container_pointer_type; + + /** + * Constructor from an accessor. + */ + Iterator (const Accessor &accessor); + + /** + * Constructor. Create the end iterator for a table. + */ + Iterator (const container_pointer_type object); + + /** + * Constructor for a particular table entry. + */ + Iterator (const container_pointer_type object, + const size_type row, + const size_type column); + + /** + * Copy constructor from a non-const iterator. + */ + Iterator(const Iterator &i); + + /** + * Constructor for an entry with a particular linear index. + */ + Iterator (const container_pointer_type container, + const std::ptrdiff_t linear_index); + }; +} + + /** * A class representing a transpose two-dimensional table, i.e. a matrix of * objects (not necessarily only numbers) in column first numbering (FORTRAN @@ -1536,6 +1734,17 @@ public: */ typedef typename AlignedVector::const_reference const_reference; + /** + * Typedef for a constant iterator that traverses the table in column-major + * order. + */ + typedef TransposeTableIterators::Iterator const_iterator; + + /** + * Typedef for an iterator that traverses the table in column-major order. + */ + typedef TransposeTableIterators::Iterator iterator; + /** * Default constructor. Set all dimensions to zero. */ @@ -1586,6 +1795,26 @@ public: */ size_type n_cols () const; + /** + * Return an iterator pointing to the first entry. + */ + iterator begin (); + + /** + * Return a constant iterator pointing to the first entry. + */ + const_iterator begin () const; + + /** + * Return an iterator pointing to one past the last entry. + */ + iterator end (); + + /** + * Return a constant iterator pointing to one past the last entry. + */ + const_iterator end () const; + protected: /** * Return a read-write reference to the element (i,j). @@ -1616,6 +1845,19 @@ protected: */ const_reference el (const size_type i, const size_type j) const; + + /** + * Make the AccessorBase class a friend so that it may directly index into + * the array. + */ + friend class TransposeTableIterators::AccessorBase; + friend class TransposeTableIterators::AccessorBase; + + /** + * Make the mutable accessor class a friend so that we can write to array + * entries with iterators. + */ + friend class TransposeTableIterators::Accessor; }; @@ -2382,6 +2624,178 @@ Table<2,T>::n_cols () const +//--------------------------------------------------------------------------- +namespace TransposeTableIterators +{ + template + inline + AccessorBase::AccessorBase() + : + container(nullptr), + linear_index(std::numeric_limits::max()) + {} + + + + template + inline + AccessorBase::AccessorBase(const container_pointer_type table) + : + container(table), + linear_index(container->values.size()) + {} + + + + template + inline + AccessorBase::AccessorBase(const AccessorBase &a) + : + container(a.container), + linear_index(a.linear_index) + {} + + + + template + inline + AccessorBase::AccessorBase(const container_pointer_type table, + const std::ptrdiff_t index) + : + container(table), + linear_index(index) + { + Assert(0 <= linear_index && linear_index < container->values.size() + 1, + ExcMessage("The current iterator points outside of the table and is " + "not the end iterator.")); + } + + + + template + inline + const T & + AccessorBase::value() const + { + Assert(0 <= linear_index && linear_index < container->values.size(), + ExcMessage("The current iterator points outside of the table.")); + return this->container->values[linear_index]; + } + + + + template + inline + typename AccessorBase::size_type + AccessorBase::AccessorBase::row() const + { + Assert(!container->empty(), ExcMessage("An empty table has no rows or columns.")); + const auto row_n = linear_index % container->n_rows(); + Assert(0 <= row_n && row_n < container->n_rows(), + ExcMessage("The current iterator points outside the table.")); + return row_n; + } + + + + template + inline + typename AccessorBase::size_type + AccessorBase::AccessorBase::column() const + { + Assert(!container->empty(), ExcMessage("An empty table has no rows or columns.")); + const auto column_n = linear_index/container->n_rows(); + Assert(0 <= column_n && column_n < container->n_cols(), + ExcMessage("The current iterator points outside the table.")); + return column_n; + } + + + template + inline + const Accessor & + Accessor::operator = (const T &t) const + { + Assert(0 <= this->linear_index && this->linear_index < this->container->values.size(), + ExcMessage("The current iterator points outside of the table.")); + this->container->values[this->linear_index] = t; + return *this; + } + + + + template + inline + const Accessor & + Accessor::operator = (T &&t) const + { + Assert(0 <= this->linear_index && this->linear_index < this->container->values.size(), + ExcMessage("The current iterator points outside of the table.")); + this->container->values[this->linear_index] = t; + return *this; + } + + + + template + inline + T & + Accessor::value() const + { + Assert(0 <= this->linear_index && this->linear_index < this->container->values.size(), + ExcMessage("The current iterator points outside of the table.")); + return this->container->values[this->linear_index]; + } + + + + template + Iterator::Iterator(const Accessor &a) + : + LinearIndexIterator, Accessor>(a) + {} + + + + template + Iterator::Iterator(const container_pointer_type table) + : + LinearIndexIterator, Accessor>(Accessor(table)) + {} + + + + template + Iterator::Iterator(const Iterator &i) + : + LinearIndexIterator, Accessor>(*i) + {} + + + + + template + Iterator::Iterator(const container_pointer_type table, + const size_type row_n, + const size_type col_n) + : + Iterator(table, table->n_rows()*col_n + row_n) + {} + + + + template + Iterator::Iterator + (const container_pointer_type table, + const std::ptrdiff_t linear_index) + : + LinearIndexIterator, Accessor> + (Accessor(table, linear_index)) + {} +} + + + //--------------------------------------------------------------------------- template inline @@ -2473,6 +2887,43 @@ TransposeTable::n_cols () const +template +inline +typename TransposeTable::iterator +TransposeTable::begin () +{ + return typename TransposeTable::iterator(this, 0, 0); +} + + + +template +inline +typename TransposeTable::const_iterator +TransposeTable::begin () const +{ + return typename TransposeTable::const_iterator(this, 0, 0); +} + + + +template +inline +typename TransposeTable::iterator +TransposeTable::end () +{ + return typename TransposeTable::iterator(this); +} + + + +template +inline +typename TransposeTable::const_iterator +TransposeTable::end () const +{ + return typename TransposeTable::const_iterator(this); +} //--------------------------------------------------------------------------- diff --git a/tests/base/transpose_table_iterators_0.cc b/tests/base/transpose_table_iterators_0.cc new file mode 100644 index 0000000000..e3ccf2f39c --- /dev/null +++ b/tests/base/transpose_table_iterators_0.cc @@ -0,0 +1,102 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2018 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, 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. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + +// Check some basic behavior of the TransposeTable iterators +#include "../tests.h" + +#include + +#include + +template +void test() +{ + Container table(3, 4); + + auto it = table.begin(); + auto t2 = it; + deallog << std::boolalpha; + deallog << boost::core::demangle(typeid(Container).name()) << std::endl; + deallog << "++it == it: " << (++table.begin() == table.begin()) << std::endl; + deallog << "++it != it: " << (++table.begin() != table.begin()) << std::endl; + deallog << "it++ == it: " << (table.begin()++ == table.begin()) << std::endl; + deallog << "it++ != it + 1: " << (table.begin()++ != table.begin() + 1) << std::endl; + deallog << "++(++it) == it + 2: " << (++(++table.begin()) == table.begin() + 2) << std::endl; + deallog << "it++ == it + 1: " << (table.begin()++ == table.begin() + 1) << std::endl; + deallog << "--(++it) == it: " << (--(++table.begin()) == table.begin()) << std::endl; + deallog << "(it--)-- == it: " << ((table.end()--)-- == table.end()) << std::endl; + deallog << "end - 5 <= end: " << ((table.end() - 5) <= table.end()) << std::endl; + deallog << "end <= end: " << (table.end() <= table.end()) << std::endl; + deallog << "end >= end: " << (table.end() >= table.end()) << std::endl; + deallog << "begin + 5 >= begin: " << ((table.begin() + 5) >= table.begin()) << std::endl; + deallog << "begin + 5 <= begin: " << ((table.begin() + 5) <= table.begin()) << std::endl; + deallog << "end - begin: " << table.end() - table.begin() << std::endl; + deallog << "begin - end: " << table.begin() - table.end() << std::endl; + deallog << "(begin + 2) - begin: " << (table.begin() + 2) - table.begin() << std::endl; + deallog << "begin + 6 < begin: " << ((table.begin() + 6) < table.begin()) << std::endl; + deallog << "begin < begin + 1: " << (table.begin() < (table.begin() + 6)) << std::endl; + deallog << "end - 5 < begin: " << ((table.end() - 6) < table.begin()) << std::endl; + t2 = it + 5; + deallog << "it+5 == (t2 = (it+5)):" << (it + 5 == t2) << std::endl; + const auto it2 = table.end() - 5; + deallog << "end - 5 position: " << it2->row() << ", " << it2->column() << std::endl; + const auto it3 = table.begin() + 5; + deallog << "begin + 5 position: " << it3->row() << ", " << it3->column() << std::endl; +} + +int main() +{ + initlog(); + + // Test an empty table + { + TransposeTable transpose_table; + Assert(transpose_table.begin() == transpose_table.end(), + ExcMessage("The beginning and end iterators should be equal for an " + "empty table.")); + TransposeTable::const_iterator begin = transpose_table.begin(); + TransposeTable::const_iterator end = transpose_table.end(); + Assert(begin == end, ExcMessage("The beginning and end const iterators should" + " be equal for an empty table.")); + } + + // test some things with accessors + { + TransposeTable table(2, 2); + + TransposeTableIterators::Accessor a3(&table, 2); + TransposeTableIterators::Accessor a4(&table, 2); + deallog << "Accessors refer to the same entry: " + << (&(TransposeTableIterators::Accessor(a4).value()) == &(a3.value())) + << std::endl; + } + + // test a non-empty rectangular table + test>(); + test>>(); + + { + // check that we can convert + TransposeTable table(3, 4); + TransposeTable::iterator it0 = table.begin(); + TransposeTable::const_iterator it1 = it0; + deallog << "converted iterators are equal: " + << (it1 == it0) + << std::endl; + } + + deallog << "OK" << std::endl; +} diff --git a/tests/base/transpose_table_iterators_0.output b/tests/base/transpose_table_iterators_0.output new file mode 100644 index 0000000000..161b855e2a --- /dev/null +++ b/tests/base/transpose_table_iterators_0.output @@ -0,0 +1,50 @@ + +DEAL::Accessors refer to the same entry: 1 +DEAL::dealii::TransposeTable +DEAL::++it == it: false +DEAL::++it != it: true +DEAL::it++ == it: true +DEAL::it++ != it + 1: true +DEAL::++(++it) == it + 2: true +DEAL::it++ == it + 1: false +DEAL::--(++it) == it: true +DEAL::(it--)-- == it: true +DEAL::end - 5 <= end: true +DEAL::end <= end: true +DEAL::end >= end: true +DEAL::begin + 5 >= begin: true +DEAL::begin + 5 <= begin: false +DEAL::end - begin: 12 +DEAL::begin - end: -12 +DEAL::(begin + 2) - begin: 2 +DEAL::begin + 6 < begin: false +DEAL::begin < begin + 1: true +DEAL::end - 5 < begin: false +DEAL::it+5 == (t2 = (it+5)):true +DEAL::end - 5 position: 1, 2 +DEAL::begin + 5 position: 2, 1 +DEAL::dealii::TransposeTable > +DEAL::++it == it: false +DEAL::++it != it: true +DEAL::it++ == it: true +DEAL::it++ != it + 1: true +DEAL::++(++it) == it + 2: true +DEAL::it++ == it + 1: false +DEAL::--(++it) == it: true +DEAL::(it--)-- == it: true +DEAL::end - 5 <= end: true +DEAL::end <= end: true +DEAL::end >= end: true +DEAL::begin + 5 >= begin: true +DEAL::begin + 5 <= begin: false +DEAL::end - begin: 12 +DEAL::begin - end: -12 +DEAL::(begin + 2) - begin: 2 +DEAL::begin + 6 < begin: false +DEAL::begin < begin + 1: true +DEAL::end - 5 < begin: false +DEAL::it+5 == (t2 = (it+5)):true +DEAL::end - 5 position: 1, 2 +DEAL::begin + 5 position: 2, 1 +DEAL::converted iterators are equal: true +DEAL::OK diff --git a/tests/base/transpose_table_iterators_1.cc b/tests/base/transpose_table_iterators_1.cc new file mode 100644 index 0000000000..25e663f98e --- /dev/null +++ b/tests/base/transpose_table_iterators_1.cc @@ -0,0 +1,95 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2018 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, 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. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + +// Check usage of TransposeTable iterators with the table itself + +#include "../tests.h" + +#include + +#include + +int main() +{ + initlog(); + + deallog << std::boolalpha; + + // test a non-empty rectangular table + TransposeTable table(3, 4); + std::iota(table.begin(), table.end(), 1.0); + for (const auto entry : table) + { + deallog << entry.row() + << ", " + << entry.column() + << ", " + << entry.value() + << std::endl; + } + + deallog << "backwards order:" << std::endl; + auto it = table.end() - 1; + for (; it >= table.begin(); --it) + { + deallog << it->row() + << ", " + << it->column() + << ", " + << it->value() + << std::endl; + } + deallog << "iterator is one before the beginning: " + << (it == table.begin() - 1) << std::endl; + + deallog << "every other entry:" << std::endl; + it = table.begin(); + for (; it < table.end(); it += 2) + { + deallog << it->row() + << ", " + << it->column() + << ", " + << it->value() + << std::endl; + } + + // print every other entry + it = table.end() - 1; + deallog << "every other entry:" << std::endl; + for (; it >= table.begin(); it -= 2) + { + deallog << it->row() + << ", " + << it->column() + << ", " + << it->value() + << std::endl; + } + + // test some type equalities + static_assert(std::is_samevalue()), double &>::value, + "The iterator value for a non-const table should not be const."); + static_assert(std::is_samevalue()), double &>::value, + "The iterator value for a non-const table should not be const."); + + const TransposeTable &ref = table; + static_assert(std::is_samevalue()), const double &>::value, + "The iterator value for a constant table should be const."); + static_assert(std::is_samevalue()), const double &>::value, + "The iterator value for a constant table should be const."); + + deallog << "OK" << std::endl; +} diff --git a/tests/base/transpose_table_iterators_1.output b/tests/base/transpose_table_iterators_1.output new file mode 100644 index 0000000000..4614620cd9 --- /dev/null +++ b/tests/base/transpose_table_iterators_1.output @@ -0,0 +1,42 @@ + +DEAL::0, 0, 1.00000 +DEAL::1, 0, 2.00000 +DEAL::2, 0, 3.00000 +DEAL::0, 1, 4.00000 +DEAL::1, 1, 5.00000 +DEAL::2, 1, 6.00000 +DEAL::0, 2, 7.00000 +DEAL::1, 2, 8.00000 +DEAL::2, 2, 9.00000 +DEAL::0, 3, 10.0000 +DEAL::1, 3, 11.0000 +DEAL::2, 3, 12.0000 +DEAL::backwards order: +DEAL::2, 3, 12.0000 +DEAL::1, 3, 11.0000 +DEAL::0, 3, 10.0000 +DEAL::2, 2, 9.00000 +DEAL::1, 2, 8.00000 +DEAL::0, 2, 7.00000 +DEAL::2, 1, 6.00000 +DEAL::1, 1, 5.00000 +DEAL::0, 1, 4.00000 +DEAL::2, 0, 3.00000 +DEAL::1, 0, 2.00000 +DEAL::0, 0, 1.00000 +DEAL::iterator is one before the beginning: true +DEAL::every other entry: +DEAL::0, 0, 1.00000 +DEAL::2, 0, 3.00000 +DEAL::1, 1, 5.00000 +DEAL::0, 2, 7.00000 +DEAL::2, 2, 9.00000 +DEAL::1, 3, 11.0000 +DEAL::every other entry: +DEAL::2, 3, 12.0000 +DEAL::0, 3, 10.0000 +DEAL::1, 2, 8.00000 +DEAL::2, 1, 6.00000 +DEAL::0, 1, 4.00000 +DEAL::1, 0, 2.00000 +DEAL::OK