Reference documentation for deal.II version 8.5.1
petsc_sparse_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2016 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 <vector>
27 
28 DEAL_II_NAMESPACE_OPEN
29 // forward declaration
30 template <typename MatrixType> class BlockMatrixBase;
31 
32 
33 namespace PETScWrappers
34 {
48  class SparseMatrix : public MatrixBase
49  {
50  public:
51 
59  struct Traits
60  {
65  static const bool zero_addition_can_be_elided = true;
66  };
67 
71  SparseMatrix ();
72 
87  SparseMatrix (const size_type m,
88  const size_type n,
89  const size_type n_nonzero_per_row,
90  const bool is_symmetric = false);
91 
108  SparseMatrix (const size_type m,
109  const size_type n,
110  const std::vector<size_type> &row_lengths,
111  const bool is_symmetric = false);
112 
130  template <typename SparsityPatternType>
131  explicit SparseMatrix (const SparsityPatternType &sparsity_pattern,
132  const bool preset_nonzero_locations = true);
133 
143  SparseMatrix &operator = (const double d);
144 
150  void reinit (const size_type m,
151  const size_type n,
152  const size_type n_nonzero_per_row,
153  const bool is_symmetric = false);
154 
160  void reinit (const size_type m,
161  const size_type n,
162  const std::vector<size_type> &row_lengths,
163  const bool is_symmetric = false);
164 
190  template <typename SparsityPatternType>
191  void reinit (const SparsityPatternType &sparsity_pattern,
192  const bool preset_nonzero_locations = true);
193 
199  virtual const MPI_Comm &get_mpi_communicator () const;
200 
215  PetscScalar matrix_norm_square (const VectorBase &v) const;
216 
225  PetscScalar matrix_scalar_product (const VectorBase &u,
226  const VectorBase &v) const;
227 
231  size_t m() const;
232 
236  size_t n() const;
237 
238  private:
239 
243  SparseMatrix(const SparseMatrix &);
248 
254  void do_reinit (const size_type m,
255  const size_type n,
256  const size_type n_nonzero_per_row,
257  const bool is_symmetric = false);
258 
262  void do_reinit (const size_type m,
263  const size_type n,
264  const std::vector<size_type> &row_lengths,
265  const bool is_symmetric = false);
266 
270  template <typename SparsityPatternType>
271  void do_reinit (const SparsityPatternType &sparsity_pattern,
272  const bool preset_nonzero_locations);
273 
278  };
279 }
280 
281 DEAL_II_NAMESPACE_CLOSE
282 
283 #endif // DEAL_II_WITH_PETSC
284 
285 /*---------------------------- petsc_sparse_matrix.h ---------------------------*/
286 
287 #endif
288 /*---------------------------- petsc_sparse_matrix.h ---------------------------*/
PetscScalar matrix_scalar_product(const VectorBase &u, const VectorBase &v) const
PetscScalar matrix_norm_square(const VectorBase &v) const
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
PetscBooleanType is_symmetric(const double tolerance=1.e-12)
void do_reinit(const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false)