16 #ifndef dealii_petsc_sparse_matrix_h 17 # define dealii_petsc_sparse_matrix_h 20 # include <deal.II/base/config.h> 22 # ifdef DEAL_II_WITH_PETSC 24 # include <deal.II/lac/exceptions.h> 25 # include <deal.II/lac/petsc_matrix_base.h> 26 # include <deal.II/lac/petsc_vector.h> 30 DEAL_II_NAMESPACE_OPEN
32 template <
typename MatrixType>
112 const std::vector<size_type> &row_lengths,
132 template <
typename SparsityPatternType>
133 explicit SparseMatrix(
const SparsityPatternType &sparsity_pattern,
134 const bool preset_nonzero_locations =
true);
178 const std::vector<size_type> &row_lengths,
206 template <
typename SparsityPatternType>
208 reinit(
const SparsityPatternType &sparsity_pattern,
209 const bool preset_nonzero_locations =
true);
216 virtual const MPI_Comm &
272 const std::vector<size_type> &row_lengths,
278 template <
typename SparsityPatternType>
280 do_reinit(
const SparsityPatternType &sparsity_pattern,
281 const bool preset_nonzero_locations);
435 const size_type n_offdiag_nonzero_per_row = 0);
466 const std::vector<size_type> &row_lengths,
468 const std::vector<size_type> &offdiag_row_lengths =
469 std::vector<size_type>());
493 template <
typename SparsityPatternType>
495 const SparsityPatternType & sparsity_pattern,
496 const std::vector<size_type> &local_rows_per_process,
497 const std::vector<size_type> &local_columns_per_process,
498 const unsigned int this_process,
499 const bool preset_nonzero_locations =
true);
538 const size_type n_offdiag_nonzero_per_row = 0);
555 const std::vector<size_type> &row_lengths,
557 const std::vector<size_type> &offdiag_row_lengths =
558 std::vector<size_type>());
579 template <
typename SparsityPatternType>
582 const SparsityPatternType & sparsity_pattern,
583 const std::vector<size_type> &local_rows_per_process,
584 const std::vector<size_type> &local_columns_per_process,
585 const unsigned int this_process,
586 const bool preset_nonzero_locations =
true);
594 template <
typename SparsityPatternType>
598 const SparsityPatternType &sparsity_pattern,
613 virtual const MPI_Comm &
626 <<
"The number of local rows " << arg1
627 <<
" must be larger than the total number of rows " 720 const size_type n_offdiag_nonzero_per_row = 0);
734 const std::vector<size_type> &row_lengths,
736 const std::vector<size_type> &offdiag_row_lengths =
737 std::vector<size_type>());
742 template <
typename SparsityPatternType>
744 do_reinit(
const SparsityPatternType & sparsity_pattern,
745 const std::vector<size_type> &local_rows_per_process,
746 const std::vector<size_type> &local_columns_per_process,
747 const unsigned int this_process,
748 const bool preset_nonzero_locations);
753 template <
typename SparsityPatternType>
757 const SparsityPatternType &sparsity_pattern);
769 inline const MPI_Comm &
777 DEAL_II_NAMESPACE_CLOSE
779 # endif // DEAL_II_WITH_PETSC static ::ExceptionBase & ExcLocalRowsTooLarge(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
void copy_from(const SparseMatrix &other)
IndexSet locally_owned_domain_indices() const
void do_reinit(const size_type m, const size_type n, const size_type local_rows, const size_type local_columns, const size_type n_nonzero_per_row, const bool is_symmetric=false, const size_type n_offdiag_nonzero_per_row=0)
virtual const MPI_Comm & get_mpi_communicator() const override
static const bool zero_addition_can_be_elided
PetscBool is_symmetric(const double tolerance=1.e-12)
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
static const bool zero_addition_can_be_elided
SparseMatrix & operator=(const double d)
void reinit(const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false)
SparseMatrix & operator=(const value_type d)
PetscScalar matrix_scalar_product(const Vector &u, const Vector &v) const
virtual const MPI_Comm & get_mpi_communicator() const override
void reinit(const MPI_Comm &communicator, const size_type m, const size_type n, const size_type local_rows, const size_type local_columns, const size_type n_nonzero_per_row, const bool is_symmetric=false, const size_type n_offdiag_nonzero_per_row=0)
IndexSet locally_owned_range_indices() const
unsigned int global_dof_index
types::global_dof_index size_type
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
PetscScalar matrix_norm_square(const Vector &v) const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
void do_reinit(const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false)