size_type
column_number(const size_type row, const size_type index) const;
+ /**
+ * Return index of column @p col in row @p row. If the column does not
+ * exist in this sparsity pattern, the returned value will be
+ * 'numbers::invalid_size_type'.
+ */
+ size_type
+ column_index(const size_type row, const size_type col) const;
+
/**
* @name Iterators
*/
// ---------------------------------------------------------------------
#include <deal.II/base/memory_consumption.h>
+#include <deal.II/base/utilities.h>
#include <deal.II/lac/dynamic_sparsity_pattern.h>
#include <deal.II/lac/sparsity_pattern.h>
}
+
+types::global_dof_index
+DynamicSparsityPattern::column_index(
+ const DynamicSparsityPattern::size_type row,
+ const DynamicSparsityPattern::size_type col) const
+{
+ Assert(row < n_rows(),
+ ExcIndexRangeType<DynamicSparsityPattern::size_type>(row,
+ 0,
+ n_rows()));
+ Assert(col < n_cols(),
+ ExcIndexRangeType<DynamicSparsityPattern::size_type>(row,
+ 0,
+ n_cols()));
+ Assert(rowset.size() == 0 || rowset.is_element(row), ExcInternalError());
+
+ const DynamicSparsityPattern::size_type local_row =
+ rowset.size() ? rowset.index_within_set(row) : row;
+
+ // now we need to do a binary search. Note that col indices are assumed to
+ // be sorted.
+ const auto &cols = lines[local_row].entries;
+ auto it = Utilities::lower_bound(cols.begin(), cols.end(), col);
+
+ if ((it != cols.end()) && (*it == col))
+ return (it - cols.begin());
+ else
+ return numbers::invalid_size_type;
+}
+
+
+
// explicit instantiations
template void
DynamicSparsityPattern::Line::add_entries(size_type *, size_type *, const bool);
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2004 - 2017 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check DynamicSparsityPattern::column_index.
+// this test is based on dynamic_sparsity_pattern_04.cc
+
+#include <deal.II/lac/dynamic_sparsity_pattern.h>
+
+#include "../tests.h"
+
+
+void
+test()
+{
+ const unsigned int N = 100;
+ DynamicSparsityPattern csp(N, N);
+ for (unsigned int i = 0; i < N; ++i)
+ for (unsigned int j = 0; j < 10; ++j)
+ csp.add(i, (i + (i + 1) * (j * j + i)) % N);
+
+ for (unsigned int i = 0; i < N; ++i)
+ for (unsigned int ind = 0; ind < csp.row_length(i); ++ind)
+ {
+ const auto j = csp.column_number(i, ind);
+ const auto val = csp.column_index(i, j);
+ AssertThrow(val == ind,
+ ExcMessage(std::to_string(val) +
+ "!=" + std::to_string(ind)));
+ }
+
+ deallog << "Ok" << std::endl;
+}
+
+
+
+int
+main()
+{
+ initlog();
+
+ test();
+ return 0;
+}