Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
petsc_block_sparse_matrix.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2004 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_petsc_block_sparse_matrix_h
16#define dealii_petsc_block_sparse_matrix_h
17
18
19#include <deal.II/base/config.h>
20
21#ifdef DEAL_II_WITH_PETSC
22
28
29# include <cmath>
30
32
33
34
35namespace PETScWrappers
36{
37 namespace MPI
38 {
66 class BlockSparseMatrix : public BlockMatrixBase<SparseMatrix>
67 {
68 public:
73
78
90
103
111 explicit BlockSparseMatrix(const Mat &);
112
116 template <size_t block_rows, size_t block_columns>
117 explicit BlockSparseMatrix(
118 const std::array<std::array<Mat, block_columns>, block_rows> &);
119
123 ~BlockSparseMatrix() override;
124
131
142 operator=(const double d);
143
157 void
158 reinit(const size_type n_block_rows, const size_type n_block_columns);
159
160
170 void
171 reinit(const std::vector<IndexSet> &rows,
172 const std::vector<IndexSet> &cols,
173 const BlockDynamicSparsityPattern &bdsp,
174 const MPI_Comm com);
175
176
180 void
181 reinit(const std::vector<IndexSet> &sizes,
182 const BlockDynamicSparsityPattern &bdsp,
183 const MPI_Comm com);
184
185
191 void
192 reinit(Mat A);
193
194
199 void
200 vmult(BlockVector &dst, const BlockVector &src) const;
201
206 void
207 vmult(BlockVector &dst, const Vector &src) const;
208
213 void
214 vmult(Vector &dst, const BlockVector &src) const;
215
220 void
221 vmult(Vector &dst, const Vector &src) const;
222
228 void
229 Tvmult(BlockVector &dst, const BlockVector &src) const;
230
235 void
236 Tvmult(BlockVector &dst, const Vector &src) const;
237
242 void
243 Tvmult(Vector &dst, const BlockVector &src) const;
244
249 void
250 Tvmult(Vector &dst, const Vector &src) const;
251
258 void
260
269 void
271
276 std::vector<IndexSet>
278
284 std::vector<IndexSet>
286
292 std::uint64_t
293 n_nonzero_elements() const;
294
299 get_mpi_communicator() const;
300
306
314 operator const Mat &() const;
315
326 Mat &
327 petsc_matrix();
328
329 private:
337
341 void
343
349 void
351 };
352
353
354
357 // ------------- inline and template functions -----------------
359 : BaseClass()
360 , petsc_nest_matrix(nullptr)
361 {}
362
363
364
367 {
368 this->reinit(A);
369 }
370
371
372
373 template <size_t block_rows, size_t block_columns>
375 const std::array<std::array<Mat, block_columns>, block_rows> &arrayA)
377 {
378 this->reinit(block_rows, block_columns);
379 this->sub_objects.reinit(block_rows, block_columns);
380 for (auto r = 0; r < block_rows; ++r)
381 for (auto c = 0; c < block_columns; ++c)
382 {
383 if (arrayA[r][c])
384 this->sub_objects[r][c] = new BlockType(arrayA[r][c]);
385 else
386 this->sub_objects[r][c] = nullptr;
387 }
388 this->collect_sizes();
389 }
390
391
392
393 inline BlockSparseMatrix &
395 {
397
398 for (size_type r = 0; r < this->n_block_rows(); ++r)
399 for (size_type c = 0; c < this->n_block_cols(); ++c)
400 this->block(r, c) = d;
401
402 return *this;
403 }
404
405
406
407 inline void
409 {
411 }
412
413
414
415 inline void
417 {
419 }
420
421
422
423 inline void
425 {
427 }
428
429
430
431 inline void
432 BlockSparseMatrix::vmult(Vector &dst, const Vector &src) const
433 {
435 }
436
437
438 inline void
440 {
442 }
443
444
445
446 inline void
448 {
450 }
451
452
453
454 inline void
456 {
458 }
459
460
461
462 inline void
464 {
466 }
467
468 } // namespace MPI
469
470} // namespace PETScWrappers
471
472
474
475
476#endif // DEAL_II_WITH_PETSC
477
478#endif // dealii_petsc_block_sparse_matrix_h
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
unsigned int n_block_rows() const
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, true > > const_iterator
void vmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
types::global_dof_index size_type
Table< 2, SmartPointer< BlockType, BlockMatrixBase< SparseMatrix > > > sub_objects
typename BlockType::value_type value_type
void vmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
unsigned int n_block_cols() const
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, false > > iterator
BlockType & block(const unsigned int row, const unsigned int column)
void vmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
void Tvmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
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 &)
void compress(VectorOperation::values operation)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
#define Assert(cond, exc)