Reference documentation for deal.II version 9.3.3
\(\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 - 2020 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
27
28# include <vector>
29
31// forward declaration
32# ifndef DOXYGEN
33template <typename MatrixType>
34class BlockMatrixBase;
35# endif
36
37namespace PETScWrappers
38{
51 class SparseMatrix : public MatrixBase
52 {
53 public:
61 struct Traits
62 {
67 static const bool zero_addition_can_be_elided = true;
68 };
69
74
90 const size_type n,
91 const size_type n_nonzero_per_row,
92 const bool is_symmetric = false);
93
111 const size_type n,
112 const std::vector<size_type> &row_lengths,
113 const bool is_symmetric = false);
114
132 template <typename SparsityPatternType>
133 explicit SparseMatrix(const SparsityPatternType &sparsity_pattern,
134 const bool preset_nonzero_locations = true);
135
146 operator=(const double d);
147
151 SparseMatrix(const SparseMatrix &) = delete;
152
157 operator=(const SparseMatrix &) = delete;
158
164 void
165 reinit(const size_type m,
166 const size_type n,
167 const size_type n_nonzero_per_row,
168 const bool is_symmetric = false);
169
175 void
176 reinit(const size_type m,
177 const size_type n,
178 const std::vector<size_type> &row_lengths,
179 const bool is_symmetric = false);
180
206 template <typename SparsityPatternType>
207 void
208 reinit(const SparsityPatternType &sparsity_pattern,
209 const bool preset_nonzero_locations = true);
210
216 virtual const MPI_Comm &
217 get_mpi_communicator() const override;
218
222 size_t
223 m() const;
224
228 size_t
229 n() const;
230
237 void
239 const SparseMatrix &B,
240 const MPI::Vector & V = MPI::Vector()) const;
241
249 void
251 const SparseMatrix &B,
252 const MPI::Vector & V = MPI::Vector()) const;
253
254 private:
260 void
261 do_reinit(const size_type m,
262 const size_type n,
263 const size_type n_nonzero_per_row,
264 const bool is_symmetric = false);
265
269 void
270 do_reinit(const size_type m,
271 const size_type n,
272 const std::vector<size_type> &row_lengths,
273 const bool is_symmetric = false);
274
278 template <typename SparsityPatternType>
279 void
280 do_reinit(const SparsityPatternType &sparsity_pattern,
281 const bool preset_nonzero_locations);
282
283 // To allow calling protected prepare_add() and prepare_set().
284 friend class BlockMatrixBase<SparseMatrix>;
285 };
286
287 namespace MPI
288 {
366 {
367 public:
372
380 struct Traits
381 {
389 static const bool zero_addition_can_be_elided = false;
390 };
391
395 SparseMatrix();
396
400 ~SparseMatrix() override;
401
426 const size_type m,
427 const size_type n,
428 const size_type local_rows,
429 const size_type local_columns,
430 const size_type n_nonzero_per_row,
431 const bool is_symmetric = false,
432 const size_type n_offdiag_nonzero_per_row = 0);
433
459 const size_type m,
460 const size_type n,
461 const size_type local_rows,
462 const size_type local_columns,
463 const std::vector<size_type> &row_lengths,
464 const bool is_symmetric = false,
465 const std::vector<size_type> &offdiag_row_lengths =
466 std::vector<size_type>());
467
490 template <typename SparsityPatternType>
492 const SparsityPatternType & sparsity_pattern,
493 const std::vector<size_type> &local_rows_per_process,
494 const std::vector<size_type> &local_columns_per_process,
495 const unsigned int this_process,
496 const bool preset_nonzero_locations = true);
497
508 operator=(const value_type d);
509
510
515 void
516 copy_from(const SparseMatrix &other);
517
527 void
529 const size_type m,
530 const size_type n,
531 const size_type local_rows,
532 const size_type local_columns,
533 const size_type n_nonzero_per_row,
534 const bool is_symmetric = false,
535 const size_type n_offdiag_nonzero_per_row = 0);
536
546 void
548 const size_type m,
549 const size_type n,
550 const size_type local_rows,
551 const size_type local_columns,
552 const std::vector<size_type> &row_lengths,
553 const bool is_symmetric = false,
554 const std::vector<size_type> &offdiag_row_lengths =
555 std::vector<size_type>());
556
576 template <typename SparsityPatternType>
577 void
579 const SparsityPatternType & sparsity_pattern,
580 const std::vector<size_type> &local_rows_per_process,
581 const std::vector<size_type> &local_columns_per_process,
582 const unsigned int this_process,
583 const bool preset_nonzero_locations = true);
584
591 template <typename SparsityPatternType>
592 void
593 reinit(const IndexSet & local_rows,
594 const IndexSet & local_columns,
595 const SparsityPatternType &sparsity_pattern,
596 const MPI_Comm & communicator);
597
603 void
604 reinit(const SparseMatrix &other);
605
610 virtual const MPI_Comm &
611 get_mpi_communicator() const override;
612
621 int,
622 int,
623 << "The number of local rows " << arg1
624 << " must be larger than the total number of rows "
625 << arg2);
627
643 PetscScalar
644 matrix_norm_square(const Vector &v) const;
645
654 PetscScalar
655 matrix_scalar_product(const Vector &u, const Vector &v) const;
656
663
671
678 void
680 const SparseMatrix &B,
681 const MPI::Vector & V = MPI::Vector()) const;
682
690 void
692 const SparseMatrix &B,
693 const MPI::Vector & V = MPI::Vector()) const;
694
695 private:
700
710 void
711 do_reinit(const size_type m,
712 const size_type n,
713 const size_type local_rows,
714 const size_type local_columns,
715 const size_type n_nonzero_per_row,
716 const bool is_symmetric = false,
717 const size_type n_offdiag_nonzero_per_row = 0);
718
726 void
727 do_reinit(const size_type m,
728 const size_type n,
729 const size_type local_rows,
730 const size_type local_columns,
731 const std::vector<size_type> &row_lengths,
732 const bool is_symmetric = false,
733 const std::vector<size_type> &offdiag_row_lengths =
734 std::vector<size_type>());
735
739 template <typename SparsityPatternType>
740 void
741 do_reinit(const SparsityPatternType & sparsity_pattern,
742 const std::vector<size_type> &local_rows_per_process,
743 const std::vector<size_type> &local_columns_per_process,
744 const unsigned int this_process,
745 const bool preset_nonzero_locations);
746
750 template <typename SparsityPatternType>
751 void
752 do_reinit(const IndexSet & local_rows,
753 const IndexSet & local_columns,
754 const SparsityPatternType &sparsity_pattern);
755
756 // To allow calling protected prepare_add() and prepare_set().
757 friend class BlockMatrixBase<SparseMatrix>;
758 };
759
760
761
762 // -------- template and inline functions ----------
763
764 inline const MPI_Comm &
766 {
767 return communicator;
768 }
769 } // namespace MPI
770} // namespace PETScWrappers
771
773
774# endif // DEAL_II_WITH_PETSC
775
776#endif
777/*--------------------------- petsc_sparse_matrix.h -------------------------*/
SparseMatrix & operator=(const value_type d)
void copy_from(const SparseMatrix &other)
virtual const MPI_Comm & get_mpi_communicator() const override
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)
PetscBool is_symmetric(const double tolerance=1.e-12)
void mmult(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)
virtual const MPI_Comm & get_mpi_communicator() const override
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
void reinit(const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false)
SparseMatrix(const SparseMatrix &)=delete
SparseMatrix & operator=(const SparseMatrix &)=delete
SparseMatrix & operator=(const double d)
#define DEAL_II_DEPRECATED
Definition: config.h:162
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
PetscScalar matrix_scalar_product(const Vector &u, const Vector &v) const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
static ::ExceptionBase & ExcLocalRowsTooLarge(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:538
PetscScalar matrix_norm_square(const Vector &v) const
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)
static const char V
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
unsigned int global_dof_index
Definition: types.h:76