From: ESeNonFossiIo Date: Wed, 13 Apr 2016 15:15:38 +0000 (+0200) Subject: PETSc ReinitHelper for BlockVector X-Git-Tag: v8.5.0-rc1~1104^2~2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=7fbc164bcccad9f1146355f92e87d7318d1ec9f8;p=dealii.git PETSc ReinitHelper for BlockVector --- diff --git a/include/deal.II/lac/block_vector.h b/include/deal.II/lac/block_vector.h index ac8a467d42..78342c74af 100644 --- a/include/deal.II/lac/block_vector.h +++ b/include/deal.II/lac/block_vector.h @@ -488,7 +488,7 @@ namespace internal template class ReinitHelper; /** - * A helper class internally used in linear_operator.h. Specialization for + * A helper class used internally in linear_operator.h. Specialization for * BlockVector. */ template diff --git a/include/deal.II/lac/petsc_block_sparse_matrix.h b/include/deal.II/lac/petsc_block_sparse_matrix.h index c4d7233448..41115858ca 100644 --- a/include/deal.II/lac/petsc_block_sparse_matrix.h +++ b/include/deal.II/lac/petsc_block_sparse_matrix.h @@ -206,6 +206,19 @@ namespace PETScWrappers void Tvmult (Vector &dst, const Vector &src) const; + /** + * Return the partitioning of the domain space of this matrix, i.e., the + * partitioning of the vectors this matrix has to be multiplied with. + */ + std::vector< size_type > locally_domain_sizes() const; + + /** + * Return the partitioning of the range space of this matrix, i.e., the + * partitioning of the vectors that are result from matrix-vector + * products. + */ + std::vector< size_type > locally_range_sizes() const; + /** * Make the clear() function in the base class visible, though it is * protected. diff --git a/include/deal.II/lac/petsc_block_vector.h b/include/deal.II/lac/petsc_block_vector.h index 9dc85b424c..e912a3367c 100644 --- a/include/deal.II/lac/petsc_block_vector.h +++ b/include/deal.II/lac/petsc_block_vector.h @@ -459,6 +459,42 @@ namespace PETScWrappers } +namespace internal +{ + namespace LinearOperator + { + template class ReinitHelper; + + /** + * A helper class used internally in linear_operator.h. Specialization for + * PETScWrappers::BlockVector. + */ + template<> + class ReinitHelper + { + public: + template + static + void reinit_range_vector (const Matrix &matrix, + PETScWrappers::BlockVector &v, + bool omit_zeroing_entries) + { + v.reinit(matrix.locally_range_sizes(), omit_zeroing_entries); + } + + template + static + void reinit_domain_vector(const Matrix &matrix, + PETScWrappers::BlockVector &v, + bool omit_zeroing_entries) + { + v.reinit(matrix.locally_domain_sizes(), omit_zeroing_entries); + } + }; + + } /* namespace LinearOperator */ +} /* namespace internal */ + DEAL_II_NAMESPACE_CLOSE #endif // DEAL_II_WITH_PETSC diff --git a/include/deal.II/lac/petsc_parallel_block_sparse_matrix.h b/include/deal.II/lac/petsc_parallel_block_sparse_matrix.h index 1f79fa60ed..fb5fb17786 100644 --- a/include/deal.II/lac/petsc_parallel_block_sparse_matrix.h +++ b/include/deal.II/lac/petsc_parallel_block_sparse_matrix.h @@ -234,6 +234,19 @@ namespace PETScWrappers */ void collect_sizes (); + /** + * Return the partitioning of the domain space of this matrix, i.e., the + * partitioning of the vectors this matrix has to be multiplied with. + */ + std::vector< IndexSet > locally_owned_domain_indices() const; + + /** + * Return the partitioning of the range space of this matrix, i.e., the + * partitioning of the vectors that are result from matrix-vector + * products. + */ + std::vector< IndexSet > locally_owned_range_indices() const; + /** * Return a reference to the MPI communicator object in use with this * matrix. diff --git a/include/deal.II/lac/petsc_parallel_block_vector.h b/include/deal.II/lac/petsc_parallel_block_vector.h index ac8fa4e1b4..9263128aab 100644 --- a/include/deal.II/lac/petsc_parallel_block_vector.h +++ b/include/deal.II/lac/petsc_parallel_block_vector.h @@ -538,6 +538,42 @@ namespace PETScWrappers } +namespace internal +{ + namespace LinearOperator + { + template class ReinitHelper; + + /** + * A helper class used internally in linear_operator.h. Specialization for + * PETScWrappers::MPI::BlockVector. + */ + template<> + class ReinitHelper + { + public: + template + static + void reinit_range_vector (const Matrix &matrix, + PETScWrappers::MPI::BlockVector &v, + bool omit_zeroing_entries) + { + v.reinit(matrix.locally_owned_range_indices(), matrix.get_mpi_communicator()); + } + + template + static + void reinit_domain_vector(const Matrix &matrix, + PETScWrappers::MPI::BlockVector &v, + bool omit_zeroing_entries) + { + v.reinit(matrix.locally_owned_domain_indices(), matrix.get_mpi_communicator()); + } + }; + + } /* namespace LinearOperator */ +} /* namespace internal */ + DEAL_II_NAMESPACE_CLOSE #endif // DEAL_II_WITH_PETSC diff --git a/include/deal.II/lac/petsc_parallel_sparse_matrix.h b/include/deal.II/lac/petsc_parallel_sparse_matrix.h index 6d17876da6..f8ced1a8b7 100644 --- a/include/deal.II/lac/petsc_parallel_sparse_matrix.h +++ b/include/deal.II/lac/petsc_parallel_sparse_matrix.h @@ -387,7 +387,7 @@ namespace PETScWrappers /** * Return the partitioning of the range space of this matrix, i.e., the - * partitioning of the vectors that are result from matrix-vector + * partitioning of the vectors that result from matrix-vector * products. */ IndexSet locally_owned_range_indices() const; diff --git a/include/deal.II/lac/petsc_parallel_vector.h b/include/deal.II/lac/petsc_parallel_vector.h index 522a45f937..3f005fea9c 100644 --- a/include/deal.II/lac/petsc_parallel_vector.h +++ b/include/deal.II/lac/petsc_parallel_vector.h @@ -589,8 +589,8 @@ namespace internal template class ReinitHelper; /** - * A helper class internally used in linear_operator.h. Specialization for - * TrilinosWrappers::MPI::Vector. + * A helper class used internally in linear_operator.h. Specialization for + * PETScWrappers::MPI::Vector. */ template<> class ReinitHelper diff --git a/include/deal.II/lac/petsc_vector.h b/include/deal.II/lac/petsc_vector.h index ba4e3c4d5d..5e113f7618 100644 --- a/include/deal.II/lac/petsc_vector.h +++ b/include/deal.II/lac/petsc_vector.h @@ -398,7 +398,7 @@ namespace internal template class ReinitHelper; /** - * A helper class internally used in linear_operator.h. Specialization for + * A helper class used internally in linear_operator.h. Specialization for * PETScWrappers::Vector. */ template<> diff --git a/include/deal.II/lac/trilinos_block_vector.h b/include/deal.II/lac/trilinos_block_vector.h index 2792388a3c..965bb14751 100644 --- a/include/deal.II/lac/trilinos_block_vector.h +++ b/include/deal.II/lac/trilinos_block_vector.h @@ -464,7 +464,7 @@ namespace internal template class ReinitHelper; /** - * A helper class internally used in linear_operator.h. Specialization for + * A helper class used internally in linear_operator.h. Specialization for * TrilinosWrappers::BlockVector. */ template<> diff --git a/include/deal.II/lac/trilinos_parallel_block_vector.h b/include/deal.II/lac/trilinos_parallel_block_vector.h index 9213f88cce..f7a9cee623 100644 --- a/include/deal.II/lac/trilinos_parallel_block_vector.h +++ b/include/deal.II/lac/trilinos_parallel_block_vector.h @@ -488,7 +488,7 @@ namespace internal template class ReinitHelper; /** - * A helper class internally used in linear_operator.h. Specialization for + * A helper class used internally in linear_operator.h. Specialization for * TrilinosWrappers::MPI::BlockVector. */ template<> diff --git a/include/deal.II/lac/trilinos_vector.h b/include/deal.II/lac/trilinos_vector.h index 990311d958..2954a48f44 100644 --- a/include/deal.II/lac/trilinos_vector.h +++ b/include/deal.II/lac/trilinos_vector.h @@ -975,7 +975,7 @@ namespace internal template class ReinitHelper; /** - * A helper class internally used in linear_operator.h. Specialization for + * A helper class used internally in linear_operator.h. Specialization for * TrilinosWrappers::MPI::Vector. */ template<> @@ -1002,7 +1002,7 @@ namespace internal }; /** - * A helper class internally used in linear_operator.h. Specialization for + * A helper class used internally in linear_operator.h. Specialization for * TrilinosWrappers::Vector. */ template<> diff --git a/source/lac/petsc_block_sparse_matrix.cc b/source/lac/petsc_block_sparse_matrix.cc index bbe6817d02..23b5be69d9 100644 --- a/source/lac/petsc_block_sparse_matrix.cc +++ b/source/lac/petsc_block_sparse_matrix.cc @@ -66,7 +66,29 @@ namespace PETScWrappers } } + std::vector + BlockSparseMatrix:: + locally_domain_sizes() const + { + std::vector< size_type > index_sets; + + for ( unsigned int i=0; in_block_cols(); ++i) + index_sets.push_back(this->block(0,i).n()); + + return index_sets; + } + std::vector + BlockSparseMatrix:: + locally_range_sizes() const + { + std::vector< size_type > index_sets; + + for ( unsigned int i=0; in_block_rows(); ++i) + index_sets.push_back(this->block(i,0).m()); + + return index_sets; + } void BlockSparseMatrix::collect_sizes () diff --git a/source/lac/petsc_parallel_block_sparse_matrix.cc b/source/lac/petsc_parallel_block_sparse_matrix.cc index 3c3555f2e0..751f9f896c 100644 --- a/source/lac/petsc_parallel_block_sparse_matrix.cc +++ b/source/lac/petsc_parallel_block_sparse_matrix.cc @@ -126,7 +126,27 @@ namespace PETScWrappers BaseClass::collect_sizes (); } + std::vector< IndexSet > + BlockSparseMatrix::locally_owned_domain_indices() const + { + std::vector< IndexSet > index_sets; + + for ( unsigned int i=0; in_block_cols(); ++i) + index_sets.push_back(this->block(0,i).locally_owned_domain_indices()); + + return index_sets; + } + std::vector< IndexSet > + BlockSparseMatrix::locally_owned_range_indices() const + { + std::vector< IndexSet > index_sets; + + for ( unsigned int i=0; in_block_rows(); ++i) + index_sets.push_back(this->block(i,0).locally_owned_range_indices()); + + return index_sets; + } const MPI_Comm & BlockSparseMatrix::get_mpi_communicator () const diff --git a/source/lac/petsc_parallel_sparse_matrix.cc b/source/lac/petsc_parallel_sparse_matrix.cc index 043a8f5e00..0888d3582a 100644 --- a/source/lac/petsc_parallel_sparse_matrix.cc +++ b/source/lac/petsc_parallel_sparse_matrix.cc @@ -899,8 +899,11 @@ namespace PETScWrappers IS *cols = nullptr; ierr = MatGetSize (matrix, &n_rows, &n_cols); + AssertThrow (ierr == 0, ExcPETScError(ierr)); ierr = MatGetOwnershipIS(matrix, rows, cols); + AssertThrow (ierr == 0, ExcPETScError(ierr)); ierr = ISGetMinMax(*rows, &min, &max); + AssertThrow (ierr == 0, ExcPETScError(ierr)); IndexSet locally_owned_domain_indices(n_rows); locally_owned_domain_indices.add_range(min, max); @@ -917,8 +920,11 @@ namespace PETScWrappers IS *cols = nullptr; ierr = MatGetSize (matrix, &n_rows, &n_cols); + AssertThrow (ierr == 0, ExcPETScError(ierr)); ierr = MatGetOwnershipIS(matrix, rows, cols); + AssertThrow (ierr == 0, ExcPETScError(ierr)); ierr = ISGetMinMax(*cols, &min, &max); + AssertThrow (ierr == 0, ExcPETScError(ierr)); IndexSet locally_owned_range_indices(n_cols); locally_owned_range_indices.add_range(min, max); diff --git a/source/lac/petsc_sparse_matrix.cc b/source/lac/petsc_sparse_matrix.cc index ec0eaa8471..750f54db0b 100644 --- a/source/lac/petsc_sparse_matrix.cc +++ b/source/lac/petsc_sparse_matrix.cc @@ -313,6 +313,7 @@ namespace PETScWrappers { PetscInt m,n; PetscErrorCode ierr = MatGetSize(matrix, &m, &n); + AssertThrow (ierr == 0, ExcPETScError(ierr)); return m; } @@ -322,6 +323,7 @@ namespace PETScWrappers { PetscInt m,n; PetscErrorCode ierr = MatGetSize(matrix, &m, &n); + AssertThrow (ierr == 0, ExcPETScError(ierr)); return n; } diff --git a/tests/lac/linear_operator_09.cc b/tests/lac/linear_operator_09.cc index 41c48cc1a7..43a3279f1b 100644 --- a/tests/lac/linear_operator_09.cc +++ b/tests/lac/linear_operator_09.cc @@ -1,6 +1,6 @@ // --------------------------------------------------------------------- // -// Copyright (C) 2015 by the deal.II authors +// Copyright (C) 2016 by the deal.II authors // // This file is part of the deal.II library. // @@ -14,7 +14,7 @@ // --------------------------------------------------------------------- // Test that it is possible to instantiate a LinearOperator object for all -// different kinds of Trilinos matrices and vectors +// different kinds of PETSc matrices and vectors // TODO: A bit more tests... #include "../tests.h" @@ -22,12 +22,21 @@ #include +// Vectors: #include #include #include #include +// Block Matrix and Vectors: +#include +#include + +#include +#include + + using namespace dealii; int main(int argc, char *argv[]) @@ -47,6 +56,15 @@ int main(int argc, char *argv[]) auto op_a = linear_operator(a); } + { + PETScWrappers::BlockSparseMatrix a; + auto op_a = linear_operator(a); + } + + { + PETScWrappers::MPI::BlockSparseMatrix a; + auto op_a = linear_operator(a); + } deallog << "OK" << std::endl;