Reference documentation for deal.II version 9.1.0-pre
petsc_block_sparse_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2017 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_block_sparse_matrix_h
17 #define dealii_petsc_block_sparse_matrix_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_PETSC
23 
24 # include <deal.II/base/table.h>
25 
26 # include <deal.II/lac/block_matrix_base.h>
27 # include <deal.II/lac/block_sparsity_pattern.h>
28 # include <deal.II/lac/exceptions.h>
29 # include <deal.II/lac/petsc_block_vector.h>
30 # include <deal.II/lac/petsc_sparse_matrix.h>
31 
32 # include <cmath>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 
37 
38 namespace PETScWrappers
39 {
40  namespace MPI
41  {
69  class BlockSparseMatrix : public BlockMatrixBase<SparseMatrix>
70  {
71  public:
76 
81 
86  using pointer = BaseClass::pointer;
87  using const_pointer = BaseClass::const_pointer;
88  using reference = BaseClass::reference;
89  using const_reference = BaseClass::const_reference;
90  using size_type = BaseClass::size_type;
93 
105  BlockSparseMatrix() = default;
106 
110  ~BlockSparseMatrix() override = default;
111 
117  operator=(const BlockSparseMatrix &);
118 
129  operator=(const double d);
130 
144  void
145  reinit(const size_type n_block_rows, const size_type n_block_columns);
146 
147 
157  void
158  reinit(const std::vector<IndexSet> & rows,
159  const std::vector<IndexSet> & cols,
160  const BlockDynamicSparsityPattern &bdsp,
161  const MPI_Comm & com);
162 
163 
167  void
168  reinit(const std::vector<IndexSet> & sizes,
169  const BlockDynamicSparsityPattern &bdsp,
170  const MPI_Comm & com);
171 
172 
173 
178  void
179  vmult(BlockVector &dst, const BlockVector &src) const;
180 
185  void
186  vmult(BlockVector &dst, const Vector &src) const;
187 
192  void
193  vmult(Vector &dst, const BlockVector &src) const;
194 
199  void
200  vmult(Vector &dst, const Vector &src) const;
201 
207  void
208  Tvmult(BlockVector &dst, const BlockVector &src) const;
209 
214  void
215  Tvmult(BlockVector &dst, const Vector &src) const;
216 
221  void
222  Tvmult(Vector &dst, const BlockVector &src) const;
223 
228  void
229  Tvmult(Vector &dst, const Vector &src) const;
230 
238  void
239  collect_sizes();
240 
245  std::vector<IndexSet>
247 
253  std::vector<IndexSet>
255 
260  const MPI_Comm &
261  get_mpi_communicator() const;
262 
268  };
269 
270 
271 
274  // ------------- inline and template functions -----------------
275 
276  inline BlockSparseMatrix &
278  {
280 
281  for (size_type r = 0; r < this->n_block_rows(); ++r)
282  for (size_type c = 0; c < this->n_block_cols(); ++c)
283  this->block(r, c) = d;
284 
285  return *this;
286  }
287 
288 
289 
290  inline void
292  {
294  }
295 
296 
297 
298  inline void
300  {
302  }
303 
304 
305 
306  inline void
308  {
310  }
311 
312 
313 
314  inline void
315  BlockSparseMatrix::vmult(Vector &dst, const Vector &src) const
316  {
318  }
319 
320 
321  inline void
323  {
325  }
326 
327 
328 
329  inline void
331  {
333  }
334 
335 
336 
337  inline void
339  {
341  }
342 
343 
344 
345  inline void
346  BlockSparseMatrix::Tvmult(Vector &dst, const Vector &src) const
347  {
349  }
350 
351  } // namespace MPI
352 
353 } // namespace PETScWrappers
354 
355 
356 DEAL_II_NAMESPACE_CLOSE
357 
358 
359 #endif // DEAL_II_WITH_PETSC
360 
361 #endif // dealii_petsc_block_sparse_matrix_h
void Tvmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
void vmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
void vmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
static::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
void Tvmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
unsigned int n_block_cols() const
void vmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
typename BlockType::value_type value_type
std::vector< IndexSet > locally_owned_range_indices() const
void Tvmult(BlockVector &dst, const BlockVector &src) const
void vmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
#define Assert(cond, exc)
Definition: exceptions.h:1227
BlockSparseMatrix & operator=(const BlockSparseMatrix &)
void reinit(const size_type n_block_rows, const size_type n_block_columns)
void Tvmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
void Tvmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
BlockType & block(const unsigned int row, const unsigned int column)
void vmult(BlockVector &dst, const BlockVector &src) const
std::vector< IndexSet > locally_owned_domain_indices() const
unsigned int n_block_rows() const