Reference documentation for deal.II version 9.0.0
petsc_parallel_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_parallel_sparse_matrix_h
17 #define dealii_petsc_parallel_sparse_matrix_h
18 
19 #include <deal.II/base/config.h>
20 
21 #ifdef DEAL_II_WITH_PETSC
22 
23 # include <deal.II/lac/exceptions.h>
24 # include <deal.II/lac/petsc_matrix_base.h>
25 # include <deal.II/lac/petsc_parallel_vector.h>
26 # include <vector>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 
31 // forward declaration
32 template <typename Matrix> class BlockMatrixBase;
33 
34 
35 namespace PETScWrappers
36 {
37  namespace MPI
38  {
39 
40 
41 
119  class SparseMatrix : public MatrixBase
120  {
121  public:
126 
134  struct Traits
135  {
143  static const bool zero_addition_can_be_elided = false;
144  };
145 
149  SparseMatrix ();
150 
154  ~SparseMatrix ();
155 
178  DEAL_II_DEPRECATED
179  SparseMatrix (const MPI_Comm &communicator,
180  const size_type m,
181  const size_type n,
182  const size_type local_rows,
183  const size_type local_columns,
184  const size_type n_nonzero_per_row,
185  const bool is_symmetric = false,
186  const size_type n_offdiag_nonzero_per_row = 0);
187 
211  DEAL_II_DEPRECATED
212  SparseMatrix (const MPI_Comm &communicator,
213  const size_type m,
214  const size_type n,
215  const size_type local_rows,
216  const size_type local_columns,
217  const std::vector<size_type> &row_lengths,
218  const bool is_symmetric = false,
219  const std::vector<size_type> &offdiag_row_lengths = std::vector<size_type>());
220 
243  template <typename SparsityPatternType>
244  SparseMatrix (const MPI_Comm &communicator,
245  const SparsityPatternType &sparsity_pattern,
246  const std::vector<size_type> &local_rows_per_process,
247  const std::vector<size_type> &local_columns_per_process,
248  const unsigned int this_process,
249  const bool preset_nonzero_locations = true);
250 
261 
262 
267  void copy_from(const SparseMatrix &other);
268 
277  DEAL_II_DEPRECATED
278  void reinit (const MPI_Comm &communicator,
279  const size_type m,
280  const size_type n,
281  const size_type local_rows,
282  const size_type local_columns,
283  const size_type n_nonzero_per_row,
284  const bool is_symmetric = false,
285  const size_type n_offdiag_nonzero_per_row = 0);
286 
295  DEAL_II_DEPRECATED
296  void reinit (const MPI_Comm &communicator,
297  const size_type m,
298  const size_type n,
299  const size_type local_rows,
300  const size_type local_columns,
301  const std::vector<size_type> &row_lengths,
302  const bool is_symmetric = false,
303  const std::vector<size_type> &offdiag_row_lengths = std::vector<size_type>());
304 
324  template <typename SparsityPatternType>
325  void reinit (const MPI_Comm &communicator,
326  const SparsityPatternType &sparsity_pattern,
327  const std::vector<size_type> &local_rows_per_process,
328  const std::vector<size_type> &local_columns_per_process,
329  const unsigned int this_process,
330  const bool preset_nonzero_locations = true);
331 
338  template <typename SparsityPatternType>
339  void reinit (const IndexSet &local_rows,
340  const IndexSet &local_columns,
341  const SparsityPatternType &sparsity_pattern,
342  const MPI_Comm &communicator);
343 
349  void reinit (const SparseMatrix &other);
350 
355  virtual const MPI_Comm &get_mpi_communicator () const;
356 
365  int, int,
366  << "The number of local rows " << arg1
367  << " must be larger than the total number of rows " << arg2);
369 
385  PetscScalar matrix_norm_square (const Vector &v) const;
386 
395  PetscScalar matrix_scalar_product (const Vector &u,
396  const Vector &v) const;
397 
403 
410 
417  void mmult (SparseMatrix &C,
418  const SparseMatrix &B,
419  const MPI::Vector &V = MPI::Vector()) const;
420 
428  void Tmmult (SparseMatrix &C,
429  const SparseMatrix &B,
430  const MPI::Vector &V = MPI::Vector()) const;
431  private:
432 
436  MPI_Comm communicator;
437 
446  DEAL_II_DEPRECATED
447  void do_reinit (const size_type m,
448  const size_type n,
449  const size_type local_rows,
450  const size_type local_columns,
451  const size_type n_nonzero_per_row,
452  const bool is_symmetric = false,
453  const size_type n_offdiag_nonzero_per_row = 0);
454 
461  DEAL_II_DEPRECATED
462  void do_reinit (const size_type m,
463  const size_type n,
464  const size_type local_rows,
465  const size_type local_columns,
466  const std::vector<size_type> &row_lengths,
467  const bool is_symmetric = false,
468  const std::vector<size_type> &offdiag_row_lengths = std::vector<size_type>());
469 
473  template <typename SparsityPatternType>
474  void do_reinit (const SparsityPatternType &sparsity_pattern,
475  const std::vector<size_type> &local_rows_per_process,
476  const std::vector<size_type> &local_columns_per_process,
477  const unsigned int this_process,
478  const bool preset_nonzero_locations);
479 
483  template <typename SparsityPatternType>
484  void do_reinit (const IndexSet &local_rows,
485  const IndexSet &local_columns,
486  const SparsityPatternType &sparsity_pattern);
487 
492  };
493 
494 
495 
496 // -------- template and inline functions ----------
497 
498  inline
499  const MPI_Comm &
501  {
502  return communicator;
503  }
504  }
505 }
506 
507 DEAL_II_NAMESPACE_CLOSE
508 
509 #endif // DEAL_II_WITH_PETSC
510 
511 /*---------------------------- petsc_parallel_sparse_matrix.h ---------------------------*/
512 
513 #endif
514 /*---------------------------- petsc_parallel_sparse_matrix.h ---------------------------*/
static ::ExceptionBase & ExcLocalRowsTooLarge(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:358
void copy_from(const SparseMatrix &other)
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
PetscBool is_symmetric(const double tolerance=1.e-12)
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
SparseMatrix & operator=(const value_type d)
PetscScalar matrix_scalar_product(const Vector &u, const Vector &v) const
unsigned int global_dof_index
Definition: types.h:88
types::global_dof_index size_type
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)
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
PetscScalar matrix_norm_square(const Vector &v) const