Add assertion to make sure that begin() and end() can only be called on a processor owning the entire matrix. Also make sure that begin() works in case the first row(s) of the matrix is(are) emtpy.
Fixes #8571.
--- /dev/null
+Fixed: Make sure that begin() and end() of PETScWrappers::MatrixBase can only be called on a processor owning all rows of the matrix. Also fix a bug which caused begin() to fail when the first row(s) of the matrix is(are) empty.
+<br>
+(Sebastian Stark, 2019/08/19)
residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const;
/**
- * Iterator starting at the first entry.
+ * Iterator starting at the first entry. This can only be called on a
+ * processor owning the entire matrix. In all other cases refer to the
+ * version of begin() taking a row number as an argument.
*/
const_iterator
begin() const;
/**
- * Final iterator.
+ * Final iterator. This can only be called on a processor owning the entire
+ * matrix. In all other cases refer to the version of end() taking a row
+ * number as an argument.
*/
const_iterator
end() const;
inline MatrixBase::const_iterator
MatrixBase::begin() const
{
- return const_iterator(this, 0, 0);
+ Assert(
+ (in_local_range(0) && in_local_range(m() - 1)),
+ ExcMessage(
+ "begin() and end() can only be called on a processor owning the entire matrix. If this is a distributed matrix, use begin(row) and end(row) instead."));
+
+ // find the first non-empty row in order to make sure that the returned
+ // iterator points to something useful
+ size_type first_nonempty_row = 0;
+ while ((first_nonempty_row < m()) && (row_length(first_nonempty_row) == 0))
+ ++first_nonempty_row;
+
+ return const_iterator(this, first_nonempty_row, 0);
}
inline MatrixBase::const_iterator
MatrixBase::end() const
{
+ Assert(
+ (in_local_range(0) && in_local_range(m() - 1)),
+ ExcMessage(
+ "begin() and end() can only be called on a processor owning the entire matrix. If this is a distributed matrix, use begin(row) and end(row) instead."));
+
return const_iterator(this, m(), 0);
}
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2004 - 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check that begin() works if the first row(s) is(are) empty
+
+#include <deal.II/lac/dynamic_sparsity_pattern.h>
+#include <deal.II/lac/petsc_sparse_matrix.h>
+
+#include "../tests.h"
+
+
+void
+test()
+{
+ // test the case that the first row is empty
+ dealii::DynamicSparsityPattern dsp_1(2);
+ dsp_1.add(1, 1);
+
+ PETScWrappers::SparseMatrix K_1(dsp_1);
+
+ if (K_1.begin()->row() == 1)
+ deallog << "OK" << std::endl;
+
+ // test the case that the entire matrix is empty
+ dealii::DynamicSparsityPattern dsp_2(2);
+
+ PETScWrappers::SparseMatrix K_2(dsp_2);
+
+ if (K_2.begin() == K_2.end())
+ deallog << "OK" << std::endl;
+}
+
+int
+main(int argc, char **argv)
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+ initlog();
+ test();
+}
--- /dev/null
+
+DEAL::OK
+DEAL::OK