Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
petsc_sparse_matrix.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_petsc_sparse_matrix_h
17 # define dealii_petsc_sparse_matrix_h
18 
19 
20 # include <deal.II/base/config.h>
21 
22 # ifdef DEAL_II_WITH_PETSC
23 
24 # include <deal.II/lac/exceptions.h>
27 
28 # include <vector>
29 
31 // forward declaration
32 # ifndef DOXYGEN
33 template <typename MatrixType>
34 class BlockMatrixBase;
35 # endif
36 
37 namespace PETScWrappers
38 {
52  class SparseMatrix : public MatrixBase
53  {
54  public:
62  struct Traits
63  {
68  static const bool zero_addition_can_be_elided = true;
69  };
70 
74  SparseMatrix();
75 
90  SparseMatrix(const size_type m,
91  const size_type n,
92  const size_type n_nonzero_per_row,
93  const bool is_symmetric = false);
94 
111  SparseMatrix(const size_type m,
112  const size_type n,
113  const std::vector<size_type> &row_lengths,
114  const bool is_symmetric = false);
115 
133  template <typename SparsityPatternType>
134  explicit SparseMatrix(const SparsityPatternType &sparsity_pattern,
135  const bool preset_nonzero_locations = true);
136 
146  SparseMatrix &
147  operator=(const double d);
148 
152  SparseMatrix(const SparseMatrix &) = delete;
153 
157  SparseMatrix &
158  operator=(const SparseMatrix &) = delete;
159 
165  void
166  reinit(const size_type m,
167  const size_type n,
168  const size_type n_nonzero_per_row,
169  const bool is_symmetric = false);
170 
176  void
177  reinit(const size_type m,
178  const size_type n,
179  const std::vector<size_type> &row_lengths,
180  const bool is_symmetric = false);
181 
207  template <typename SparsityPatternType>
208  void
209  reinit(const SparsityPatternType &sparsity_pattern,
210  const bool preset_nonzero_locations = true);
211 
217  virtual const MPI_Comm &
218  get_mpi_communicator() const override;
219 
223  size_t
224  m() const;
225 
229  size_t
230  n() const;
231 
238  void
239  mmult(SparseMatrix & C,
240  const SparseMatrix &B,
241  const MPI::Vector & V = MPI::Vector()) const;
242 
250  void
252  const SparseMatrix &B,
253  const MPI::Vector & V = MPI::Vector()) const;
254 
255  private:
261  void
262  do_reinit(const size_type m,
263  const size_type n,
264  const size_type n_nonzero_per_row,
265  const bool is_symmetric = false);
266 
270  void
271  do_reinit(const size_type m,
272  const size_type n,
273  const std::vector<size_type> &row_lengths,
274  const bool is_symmetric = false);
275 
279  template <typename SparsityPatternType>
280  void
281  do_reinit(const SparsityPatternType &sparsity_pattern,
282  const bool preset_nonzero_locations);
283 
284  // To allow calling protected prepare_add() and prepare_set().
286  };
287 
288  namespace MPI
289  {
367  class SparseMatrix : public MatrixBase
368  {
369  public:
374 
382  struct Traits
383  {
391  static const bool zero_addition_can_be_elided = false;
392  };
393 
397  SparseMatrix();
398 
402  ~SparseMatrix() override;
403 
428  const size_type m,
429  const size_type n,
430  const size_type local_rows,
431  const size_type local_columns,
432  const size_type n_nonzero_per_row,
433  const bool is_symmetric = false,
434  const size_type n_offdiag_nonzero_per_row = 0);
435 
461  const size_type m,
462  const size_type n,
463  const size_type local_rows,
464  const size_type local_columns,
465  const std::vector<size_type> &row_lengths,
466  const bool is_symmetric = false,
467  const std::vector<size_type> &offdiag_row_lengths =
468  std::vector<size_type>());
469 
492  template <typename SparsityPatternType>
494  const SparsityPatternType & sparsity_pattern,
495  const std::vector<size_type> &local_rows_per_process,
496  const std::vector<size_type> &local_columns_per_process,
497  const unsigned int this_process,
498  const bool preset_nonzero_locations = true);
499 
509  SparseMatrix &
510  operator=(const value_type d);
511 
512 
517  void
518  copy_from(const SparseMatrix &other);
519 
529  void
531  const size_type m,
532  const size_type n,
533  const size_type local_rows,
534  const size_type local_columns,
535  const size_type n_nonzero_per_row,
536  const bool is_symmetric = false,
537  const size_type n_offdiag_nonzero_per_row = 0);
538 
548  void
549  reinit(const MPI_Comm & communicator,
550  const size_type m,
551  const size_type n,
552  const size_type local_rows,
553  const size_type local_columns,
554  const std::vector<size_type> &row_lengths,
555  const bool is_symmetric = false,
556  const std::vector<size_type> &offdiag_row_lengths =
557  std::vector<size_type>());
558 
578  template <typename SparsityPatternType>
579  void
580  reinit(const MPI_Comm & communicator,
581  const SparsityPatternType & sparsity_pattern,
582  const std::vector<size_type> &local_rows_per_process,
583  const std::vector<size_type> &local_columns_per_process,
584  const unsigned int this_process,
585  const bool preset_nonzero_locations = true);
586 
593  template <typename SparsityPatternType>
594  void
595  reinit(const IndexSet & local_rows,
596  const IndexSet & local_columns,
597  const SparsityPatternType &sparsity_pattern,
598  const MPI_Comm & communicator);
599 
605  void
606  reinit(const SparseMatrix &other);
607 
612  virtual const MPI_Comm &
613  get_mpi_communicator() const override;
614 
623  int,
624  int,
625  << "The number of local rows " << arg1
626  << " must be larger than the total number of rows "
627  << arg2);
629 
645  PetscScalar
646  matrix_norm_square(const Vector &v) const;
647 
656  PetscScalar
657  matrix_scalar_product(const Vector &u, const Vector &v) const;
658 
663  IndexSet
665 
671  IndexSet
673 
680  void
681  mmult(SparseMatrix & C,
682  const SparseMatrix &B,
683  const MPI::Vector & V = MPI::Vector()) const;
684 
692  void
694  const SparseMatrix &B,
695  const MPI::Vector & V = MPI::Vector()) const;
696 
697  private:
702 
712  void
713  do_reinit(const size_type m,
714  const size_type n,
715  const size_type local_rows,
716  const size_type local_columns,
717  const size_type n_nonzero_per_row,
718  const bool is_symmetric = false,
719  const size_type n_offdiag_nonzero_per_row = 0);
720 
728  void
729  do_reinit(const size_type m,
730  const size_type n,
731  const size_type local_rows,
732  const size_type local_columns,
733  const std::vector<size_type> &row_lengths,
734  const bool is_symmetric = false,
735  const std::vector<size_type> &offdiag_row_lengths =
736  std::vector<size_type>());
737 
741  template <typename SparsityPatternType>
742  void
743  do_reinit(const SparsityPatternType & sparsity_pattern,
744  const std::vector<size_type> &local_rows_per_process,
745  const std::vector<size_type> &local_columns_per_process,
746  const unsigned int this_process,
747  const bool preset_nonzero_locations);
748 
752  template <typename SparsityPatternType>
753  void
754  do_reinit(const IndexSet & local_rows,
755  const IndexSet & local_columns,
756  const SparsityPatternType &sparsity_pattern);
757 
758  // To allow calling protected prepare_add() and prepare_set().
760  };
761 
762 
763 
764  // -------- template and inline functions ----------
765 
766  inline const MPI_Comm &
768  {
769  return communicator;
770  }
771  } // namespace MPI
772 } // namespace PETScWrappers
773 
775 
776 # endif // DEAL_II_WITH_PETSC
777 
778 #endif
779 /*--------------------------- petsc_sparse_matrix.h -------------------------*/
IndexSet
Definition: index_set.h:74
PETScWrappers::SparseMatrix
Definition: petsc_sparse_matrix.h:52
PETScWrappers::MPI::SparseMatrix::Traits::zero_addition_can_be_elided
static const bool zero_addition_can_be_elided
Definition: petsc_sparse_matrix.h:391
PETScWrappers::MPI::SparseMatrix::get_mpi_communicator
virtual const MPI_Comm & get_mpi_communicator() const override
Definition: petsc_sparse_matrix.h:767
PETScWrappers::SparseMatrix::do_reinit
void do_reinit(const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false)
Definition: petsc_sparse_matrix.cc:138
PETScWrappers
Definition: petsc_block_sparse_matrix.h:36
PETScWrappers::SparseMatrix::get_mpi_communicator
virtual const MPI_Comm & get_mpi_communicator() const override
Definition: petsc_sparse_matrix.cc:126
PETScWrappers::MPI::SparseMatrix
Definition: petsc_sparse_matrix.h:367
PETScWrappers::MPI::SparseMatrix::operator=
SparseMatrix & operator=(const value_type d)
Definition: petsc_parallel_sparse_matrix.cc:130
LAPACKSupport::V
static const char V
Definition: lapack_support.h:175
PETScWrappers::SparseMatrix::Traits
Definition: petsc_sparse_matrix.h:62
PETScWrappers::MPI::SparseMatrix::Traits
Definition: petsc_sparse_matrix.h:382
PETScWrappers::MatrixBase::n
size_type n() const
Definition: petsc_matrix_base.cc:250
petsc_matrix_base.h
PETScWrappers::SparseMatrix::operator=
SparseMatrix & operator=(const double d)
Definition: petsc_sparse_matrix.cc:70
PETScWrappers::MPI::SparseMatrix::do_reinit
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)
Definition: petsc_parallel_sparse_matrix.cc:245
PETScWrappers::MPI::SparseMatrix::ExcLocalRowsTooLarge
static ::ExceptionBase & ExcLocalRowsTooLarge(int arg1, int arg2)
Physics::Elasticity::Kinematics::C
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
MPI_Comm
exceptions.h
PETScWrappers::MatrixBase::is_symmetric
PetscBool is_symmetric(const double tolerance=1.e-12)
Definition: petsc_matrix_base.cc:620
PETScWrappers::MPI::SparseMatrix::SparseMatrix
SparseMatrix()
Definition: petsc_parallel_sparse_matrix.cc:34
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
PETScWrappers::MPI::SparseMatrix::mmult
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
Definition: petsc_parallel_sparse_matrix.cc:749
PETScWrappers::SparseMatrix::SparseMatrix
SparseMatrix()
Definition: petsc_sparse_matrix.cc:30
PETScWrappers::MPI::SparseMatrix::Tmmult
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
Definition: petsc_parallel_sparse_matrix.cc:760
PETScWrappers::MPI::SparseMatrix::matrix_scalar_product
PetscScalar matrix_scalar_product(const Vector &u, const Vector &v) const
Definition: petsc_parallel_sparse_matrix.cc:688
PETScWrappers::MPI::SparseMatrix::~SparseMatrix
~SparseMatrix() override
Definition: petsc_parallel_sparse_matrix.cc:47
BlockMatrixBase
Definition: affine_constraints.h:1903
types::global_dof_index
unsigned int global_dof_index
Definition: types.h:76
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
PETScWrappers::MatrixBase
Definition: petsc_matrix_base.h:284
DEAL_II_DEPRECATED
#define DEAL_II_DEPRECATED
Definition: config.h:98
PETScWrappers::SparseMatrix::reinit
void reinit(const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false)
Definition: petsc_sparse_matrix.cc:79
PETScWrappers::SparseMatrix::mmult
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
Definition: petsc_sparse_matrix.cc:268
PETScWrappers::MPI::SparseMatrix::copy_from
void copy_from(const SparseMatrix &other)
Definition: petsc_parallel_sparse_matrix.cc:137
PETScWrappers::MPI::Vector
Definition: petsc_vector.h:158
PETScWrappers::SparseMatrix::Traits::zero_addition_can_be_elided
static const bool zero_addition_can_be_elided
Definition: petsc_sparse_matrix.h:68
PETScWrappers::MPI::SparseMatrix::locally_owned_range_indices
IndexSet locally_owned_range_indices() const
Definition: petsc_parallel_sparse_matrix.cc:723
unsigned int
PETScWrappers::SparseMatrix::n
size_t n() const
Definition: petsc_sparse_matrix.cc:258
PETScWrappers::MPI::SparseMatrix::size_type
types::global_dof_index size_type
Definition: petsc_sparse_matrix.h:373
PETScWrappers::MatrixBase::value_type
PetscScalar value_type
Definition: petsc_matrix_base.h:300
PETScWrappers::SparseMatrix::m
size_t m() const
Definition: petsc_sparse_matrix.cc:248
PETScWrappers::MPI::SparseMatrix::locally_owned_domain_indices
IndexSet locally_owned_domain_indices() const
Definition: petsc_parallel_sparse_matrix.cc:697
PETScWrappers::MPI::SparseMatrix::communicator
MPI_Comm communicator
Definition: petsc_sparse_matrix.h:701
PETScWrappers::MatrixBase::m
size_type m() const
Definition: petsc_matrix_base.cc:237
petsc_vector.h
config.h
PETScWrappers::SparseMatrix::Tmmult
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
Definition: petsc_sparse_matrix.cc:279
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
PETScWrappers::MPI::SparseMatrix::matrix_norm_square
PetscScalar matrix_norm_square(const Vector &v) const
Definition: petsc_parallel_sparse_matrix.cc:679
DeclException2
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:541
PETScWrappers::MPI::SparseMatrix::reinit
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)
Definition: petsc_parallel_sparse_matrix.cc:150