--- /dev/null
+New: Add DynamicSparsityPattern::nonempty_cols()/DynamicSparsityPattern::nonempty_rows() to return columns/rows
+stored in the sparsity pattern.
+<br>
+(Denis Davydov, 2019/03/20)
const IndexSet &
row_index_set() const;
+ /**
+ * Return the IndexSet that contains entries for all columns in which at least
+ * one element exists in this sparsity pattern.
+ *
+ * @note In a parallel context, this only considers the locally stored rows.
+ */
+ IndexSet
+ nonempty_cols() const;
+
+ /**
+ * Return the IndexSet that contains entries for all rows in which at least
+ * one element exists in this sparsity pattern.
+ *
+ * @note In a parallel context, this only considers the locally stored rows.
+ */
+ IndexSet
+ nonempty_rows() const;
+
/**
* return whether this object stores only those entries that have been added
* explicitly, or if the sparsity pattern contains elements that have been
#include <cmath>
#include <functional>
#include <numeric>
+#include <set>
DEAL_II_NAMESPACE_OPEN
}
+
+IndexSet
+DynamicSparsityPattern::nonempty_cols() const
+{
+ std::set<types::global_dof_index> cols;
+ for (const auto &line : lines)
+ cols.insert(line.entries.begin(), line.entries.end());
+
+ IndexSet res(this->n_cols());
+ res.add_indices(cols.begin(), cols.end());
+ return res;
+}
+
+
+
+IndexSet
+DynamicSparsityPattern::nonempty_rows() const
+{
+ const IndexSet all_rows = complete_index_set(this->n_rows());
+ const IndexSet &locally_stored_rows = rowset.size() == 0 ? all_rows : rowset;
+
+ std::vector<types::global_dof_index> rows;
+ auto line = lines.begin();
+ AssertDimension(locally_stored_rows.n_elements(), lines.size());
+ for (const auto &row : locally_stored_rows)
+ {
+ if (line->entries.size() > 0)
+ rows.push_back(row);
+
+ ++line;
+ }
+
+ IndexSet res(this->n_rows());
+ res.add_indices(rows.begin(), rows.end());
+ return res;
+}
+
+
+
DynamicSparsityPattern::size_type
DynamicSparsityPattern::memory_consumption() const
{
--- /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::nonempty_columns()
+// and nonempty_rows()
+
+#include <deal.II/lac/dynamic_sparsity_pattern.h>
+
+#include "../tests.h"
+
+
+void
+test()
+{
+ const unsigned int N = 100;
+ DynamicSparsityPattern dsp(N, 2 * N);
+
+ for (unsigned int i = 0; i < 5; ++i)
+ for (unsigned int j = (i == 0 ? 0 : i - 1); j < i + 1; ++j)
+ dsp.add(i, j);
+
+ dsp.add(50, 50);
+ dsp.add(51, 51);
+ dsp.add(90, 52);
+
+ dsp.add(N - 1, 2 * N - 1);
+
+ const IndexSet cols = dsp.nonempty_cols();
+ deallog << "Columns:" << std::endl;
+ cols.print(deallog.get_file_stream());
+
+ const IndexSet rows = dsp.nonempty_rows();
+ deallog << "Rows:" << std::endl;
+ rows.print(deallog.get_file_stream());
+}
+
+
+
+int
+main()
+{
+ initlog();
+
+ test();
+ return 0;
+}
--- /dev/null
+
+DEAL::Columns:
+{[0,4], [50,52], 199}
+DEAL::Rows:
+{[0,4], [50,51], 90, 99}