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_block_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_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
29
30# include <cmath>
31
33
34
35
37{
38 namespace MPI
39 {
66 class BlockSparseMatrix : public BlockMatrixBase<SparseMatrix>
67 {
68 public:
73
78
90
102 BlockSparseMatrix() = default;
103
107 ~BlockSparseMatrix() override = default;
108
115
126 operator=(const double d);
127
141 void
142 reinit(const size_type n_block_rows, const size_type n_block_columns);
143
144
154 void
155 reinit(const std::vector<IndexSet> & rows,
156 const std::vector<IndexSet> & cols,
157 const BlockDynamicSparsityPattern &bdsp,
158 const MPI_Comm & com);
159
160
164 void
165 reinit(const std::vector<IndexSet> & sizes,
166 const BlockDynamicSparsityPattern &bdsp,
167 const MPI_Comm & com);
168
169
170
175 void
176 vmult(BlockVector &dst, const BlockVector &src) const;
177
182 void
183 vmult(BlockVector &dst, const Vector &src) const;
184
189 void
190 vmult(Vector &dst, const BlockVector &src) const;
191
196 void
197 vmult(Vector &dst, const Vector &src) const;
198
204 void
205 Tvmult(BlockVector &dst, const BlockVector &src) const;
206
211 void
212 Tvmult(BlockVector &dst, const Vector &src) const;
213
218 void
219 Tvmult(Vector &dst, const BlockVector &src) const;
220
225 void
226 Tvmult(Vector &dst, const Vector &src) const;
227
235 void
237
242 std::vector<IndexSet>
244
250 std::vector<IndexSet>
252
257 const MPI_Comm &
258 get_mpi_communicator() const;
259
265 };
266
267
268
271 // ------------- inline and template functions -----------------
272
273 inline BlockSparseMatrix &
275 {
277
278 for (size_type r = 0; r < this->n_block_rows(); ++r)
279 for (size_type c = 0; c < this->n_block_cols(); ++c)
280 this->block(r, c) = d;
281
282 return *this;
283 }
284
285
286
287 inline void
289 {
291 }
292
293
294
295 inline void
297 {
299 }
300
301
302
303 inline void
305 {
307 }
308
309
310
311 inline void
312 BlockSparseMatrix::vmult(Vector &dst, const Vector &src) const
313 {
315 }
316
317
318 inline void
320 {
322 }
323
324
325
326 inline void
328 {
330 }
331
332
333
334 inline void
336 {
338 }
339
340
341
342 inline void
344 {
346 }
347
348 } // namespace MPI
349
350} // namespace PETScWrappers
351
352
354
355
356#endif // DEAL_II_WITH_PETSC
357
358#endif // dealii_petsc_block_sparse_matrix_h
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
void Tvmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
void vmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
void Tvmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
void Tvmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
void vmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
void vmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
#define Assert(cond, exc)
Definition: exceptions.h:1465
void vmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
void Tvmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
unsigned int n_block_rows() const
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, true > > const_iterator
types::global_dof_index size_type
typename BlockType::value_type value_type
unsigned int n_block_cols() const
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, false > > iterator
BlockType & block(const unsigned int row, const unsigned int column)
void Tvmult(BlockVector &dst, const BlockVector &src) const
void vmult(BlockVector &dst, const BlockVector &src) const
void reinit(const size_type n_block_rows, const size_type n_block_columns)
BlockSparseMatrix & operator=(const BlockSparseMatrix &)
std::vector< IndexSet > locally_owned_domain_indices() const
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)