Reference documentation for deal.II version 9.0.0
petsc_sparse_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2018 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 at
12 // the top level of the deal.II distribution.
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>
25 # include <deal.II/lac/petsc_matrix_base.h>
26 # include <deal.II/lac/petsc_parallel_vector.h>
27 # include <vector>
28 
29 DEAL_II_NAMESPACE_OPEN
30 // forward declaration
31 template <typename MatrixType> class BlockMatrixBase;
32 
33 
34 namespace PETScWrappers
35 {
49  class SparseMatrix : public MatrixBase
50  {
51  public:
52 
60  struct Traits
61  {
66  static const bool zero_addition_can_be_elided = true;
67  };
68 
72  SparseMatrix ();
73 
88  SparseMatrix (const size_type m,
89  const size_type n,
90  const size_type n_nonzero_per_row,
91  const bool is_symmetric = false);
92 
109  SparseMatrix (const size_type m,
110  const size_type n,
111  const std::vector<size_type> &row_lengths,
112  const bool is_symmetric = false);
113 
131  template <typename SparsityPatternType>
132  explicit SparseMatrix (const SparsityPatternType &sparsity_pattern,
133  const bool preset_nonzero_locations = true);
134 
144  SparseMatrix &operator = (const double d);
145 
151  void reinit (const size_type m,
152  const size_type n,
153  const size_type n_nonzero_per_row,
154  const bool is_symmetric = false);
155 
161  void reinit (const size_type m,
162  const size_type n,
163  const std::vector<size_type> &row_lengths,
164  const bool is_symmetric = false);
165 
191  template <typename SparsityPatternType>
192  void reinit (const SparsityPatternType &sparsity_pattern,
193  const bool preset_nonzero_locations = true);
194 
200  virtual const MPI_Comm &get_mpi_communicator () const;
201 
205  size_t m() const;
206 
210  size_t n() const;
211 
218  void mmult (SparseMatrix &C,
219  const SparseMatrix &B,
220  const MPI::Vector &V = MPI::Vector()) const;
221 
229  void Tmmult (SparseMatrix &C,
230  const SparseMatrix &B,
231  const MPI::Vector &V = MPI::Vector()) const;
232  private:
233 
237  SparseMatrix(const SparseMatrix &) = delete;
241  SparseMatrix &operator= (const SparseMatrix &) = delete;
242 
248  void do_reinit (const size_type m,
249  const size_type n,
250  const size_type n_nonzero_per_row,
251  const bool is_symmetric = false);
252 
256  void do_reinit (const size_type m,
257  const size_type n,
258  const std::vector<size_type> &row_lengths,
259  const bool is_symmetric = false);
260 
264  template <typename SparsityPatternType>
265  void do_reinit (const SparsityPatternType &sparsity_pattern,
266  const bool preset_nonzero_locations);
267 
272  };
273 }
274 
275 DEAL_II_NAMESPACE_CLOSE
276 
277 #endif // DEAL_II_WITH_PETSC
278 
279 /*---------------------------- petsc_sparse_matrix.h ---------------------------*/
280 
281 #endif
282 /*---------------------------- petsc_sparse_matrix.h ---------------------------*/
PetscBool is_symmetric(const double tolerance=1.e-12)
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)
virtual const MPI_Comm & get_mpi_communicator() const
types::global_dof_index size_type
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) 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)