Reference documentation for deal.II version 8.5.1
petsc_parallel_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_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 
175  SparseMatrix (const MPI_Comm &communicator,
176  const size_type m,
177  const size_type n,
178  const size_type local_rows,
179  const size_type local_columns,
180  const size_type n_nonzero_per_row,
181  const bool is_symmetric = false,
182  const size_type n_offdiag_nonzero_per_row = 0);
183 
204  SparseMatrix (const MPI_Comm &communicator,
205  const size_type m,
206  const size_type n,
207  const size_type local_rows,
208  const size_type local_columns,
209  const std::vector<size_type> &row_lengths,
210  const bool is_symmetric = false,
211  const std::vector<size_type> &offdiag_row_lengths = std::vector<size_type>());
212 
235  template <typename SparsityPatternType>
236  SparseMatrix (const MPI_Comm &communicator,
237  const SparsityPatternType &sparsity_pattern,
238  const std::vector<size_type> &local_rows_per_process,
239  const std::vector<size_type> &local_columns_per_process,
240  const unsigned int this_process,
241  const bool preset_nonzero_locations = true);
242 
253 
254 
259  void copy_from(const SparseMatrix &other);
260 
266  void reinit (const MPI_Comm &communicator,
267  const size_type m,
268  const size_type n,
269  const size_type local_rows,
270  const size_type local_columns,
271  const size_type n_nonzero_per_row,
272  const bool is_symmetric = false,
273  const size_type n_offdiag_nonzero_per_row = 0);
274 
280  void reinit (const MPI_Comm &communicator,
281  const size_type m,
282  const size_type n,
283  const size_type local_rows,
284  const size_type local_columns,
285  const std::vector<size_type> &row_lengths,
286  const bool is_symmetric = false,
287  const std::vector<size_type> &offdiag_row_lengths = std::vector<size_type>());
288 
308  template <typename SparsityPatternType>
309  void reinit (const MPI_Comm &communicator,
310  const SparsityPatternType &sparsity_pattern,
311  const std::vector<size_type> &local_rows_per_process,
312  const std::vector<size_type> &local_columns_per_process,
313  const unsigned int this_process,
314  const bool preset_nonzero_locations = true);
315 
322  template <typename SparsityPatternType>
323  void reinit (const IndexSet &local_rows,
324  const IndexSet &local_columns,
325  const SparsityPatternType &sparsity_pattern,
326  const MPI_Comm &communicator);
327 
333  void reinit (const SparseMatrix &other);
334 
339  virtual const MPI_Comm &get_mpi_communicator () const;
340 
349  int, int,
350  << "The number of local rows " << arg1
351  << " must be larger than the total number of rows " << arg2);
353 
369  PetscScalar matrix_norm_square (const Vector &v) const;
370 
379  PetscScalar matrix_scalar_product (const Vector &u,
380  const Vector &v) const;
381 
387 
394 
395  private:
396 
400  MPI_Comm communicator;
401 
407  void do_reinit (const size_type m,
408  const size_type n,
409  const size_type local_rows,
410  const size_type local_columns,
411  const size_type n_nonzero_per_row,
412  const bool is_symmetric = false,
413  const size_type n_offdiag_nonzero_per_row = 0);
414 
418  void do_reinit (const size_type m,
419  const size_type n,
420  const size_type local_rows,
421  const size_type local_columns,
422  const std::vector<size_type> &row_lengths,
423  const bool is_symmetric = false,
424  const std::vector<size_type> &offdiag_row_lengths = std::vector<size_type>());
425 
429  template <typename SparsityPatternType>
430  void do_reinit (const SparsityPatternType &sparsity_pattern,
431  const std::vector<size_type> &local_rows_per_process,
432  const std::vector<size_type> &local_columns_per_process,
433  const unsigned int this_process,
434  const bool preset_nonzero_locations);
435 
439  template <typename SparsityPatternType>
440  void do_reinit (const IndexSet &local_rows,
441  const IndexSet &local_columns,
442  const SparsityPatternType &sparsity_pattern);
443 
448  };
449 
450 
451 
452 // -------- template and inline functions ----------
453 
454  inline
455  const MPI_Comm &
457  {
458  return communicator;
459  }
460  }
461 }
462 
463 DEAL_II_NAMESPACE_CLOSE
464 
465 #endif // DEAL_II_WITH_PETSC
466 
467 /*---------------------------- petsc_parallel_sparse_matrix.h ---------------------------*/
468 
469 #endif
470 /*---------------------------- petsc_parallel_sparse_matrix.h ---------------------------*/
static ::ExceptionBase & ExcLocalRowsTooLarge(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
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
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)
PetscBooleanType is_symmetric(const double tolerance=1.e-12)
PetscScalar matrix_norm_square(const Vector &v) const