Reference documentation for deal.II version 9.4.1
\(\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
block_sparse_matrix.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2000 - 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_block_sparse_matrix_h
17#define dealii_block_sparse_matrix_h
18
19
20#include <deal.II/base/config.h>
21
27
28#include <cmath>
29
31
32
48template <typename number>
49class BlockSparseMatrix : public BlockMatrixBase<SparseMatrix<number>>
50{
51public:
56
61
66 using pointer = typename BaseClass::pointer;
71 using iterator = typename BaseClass::iterator;
73
89 BlockSparseMatrix() = default;
90
104
108 virtual ~BlockSparseMatrix() override;
109
110
111
118
128 operator=(const double d);
129
138 void
140
155 virtual void
156 reinit(const BlockSparsityPattern &sparsity);
158
167 bool
168 empty() const;
169
174 get_row_length(const size_type row) const;
175
183
190 n_actually_nonzero_elements(const double threshold = 0.0) const;
191
202
207 std::size_t
210
219 template <typename block_number>
220 void
222 const BlockVector<block_number> &src) const;
223
228 template <typename block_number, typename nonblock_number>
229 void
231 const Vector<nonblock_number> &src) const;
232
237 template <typename block_number, typename nonblock_number>
238 void
240 const BlockVector<block_number> &src) const;
241
246 template <typename nonblock_number>
247 void
249
255 template <typename block_number>
256 void
258 const BlockVector<block_number> &src) const;
259
264 template <typename block_number, typename nonblock_number>
265 void
267 const Vector<nonblock_number> &src) const;
268
273 template <typename block_number, typename nonblock_number>
274 void
276 const BlockVector<block_number> &src) const;
277
282 template <typename nonblock_number>
283 void
285 const Vector<nonblock_number> &src) const;
287
299 template <class BlockVectorType>
300 void
301 precondition_Jacobi(BlockVectorType & dst,
302 const BlockVectorType &src,
303 const number omega = 1.) const;
304
310 template <typename number2>
311 void
313 const Vector<number2> &src,
314 const number omega = 1.) const;
316
341 void
342 print_formatted(std::ostream & out,
343 const unsigned int precision = 3,
344 const bool scientific = true,
345 const unsigned int width = 0,
346 const char * zero_string = " ",
347 const double denominator = 1.) const;
349
359
360private:
368};
369
370
371
373/* ------------------------- Template functions ---------------------- */
374
375
376
377template <typename number>
380{
382
383 for (size_type r = 0; r < this->n_block_rows(); ++r)
384 for (size_type c = 0; c < this->n_block_cols(); ++c)
385 this->block(r, c) = d;
386
387 return *this;
388}
389
390
391
392template <typename number>
393template <typename block_number>
394inline void
396 const BlockVector<block_number> &src) const
397{
398 BaseClass::vmult_block_block(dst, src);
399}
400
401
402
403template <typename number>
404template <typename block_number, typename nonblock_number>
405inline void
407 const Vector<nonblock_number> &src) const
408{
409 BaseClass::vmult_block_nonblock(dst, src);
410}
411
412
413
414template <typename number>
415template <typename block_number, typename nonblock_number>
416inline void
418 const BlockVector<block_number> &src) const
419{
420 BaseClass::vmult_nonblock_block(dst, src);
421}
422
423
424
425template <typename number>
426template <typename nonblock_number>
427inline void
429 const Vector<nonblock_number> &src) const
430{
431 BaseClass::vmult_nonblock_nonblock(dst, src);
432}
433
434
435
436template <typename number>
437template <typename block_number>
438inline void
440 const BlockVector<block_number> &src) const
441{
442 BaseClass::Tvmult_block_block(dst, src);
443}
444
445
446
447template <typename number>
448template <typename block_number, typename nonblock_number>
449inline void
451 const Vector<nonblock_number> &src) const
452{
453 BaseClass::Tvmult_block_nonblock(dst, src);
454}
455
456
457
458template <typename number>
459template <typename block_number, typename nonblock_number>
460inline void
462 const BlockVector<block_number> &src) const
463{
464 BaseClass::Tvmult_nonblock_block(dst, src);
465}
466
467
468
469template <typename number>
470template <typename nonblock_number>
471inline void
473 const Vector<nonblock_number> &src) const
474{
475 BaseClass::Tvmult_nonblock_nonblock(dst, src);
476}
477
478
479
480template <typename number>
481template <class BlockVectorType>
482inline void
484 const BlockVectorType &src,
485 const number omega) const
486{
487 Assert(this->n_block_rows() == this->n_block_cols(), ExcNotQuadratic());
488 Assert(dst.n_blocks() == this->n_block_rows(),
489 ExcDimensionMismatch(dst.n_blocks(), this->n_block_rows()));
490 Assert(src.n_blocks() == this->n_block_cols(),
491 ExcDimensionMismatch(src.n_blocks(), this->n_block_cols()));
492
493 // do a diagonal preconditioning. uses only
494 // the diagonal blocks of the matrix
495 for (size_type i = 0; i < this->n_block_rows(); ++i)
496 this->block(i, i).precondition_Jacobi(dst.block(i), src.block(i), omega);
497}
498
499
500
501template <typename number>
502template <typename number2>
503inline void
505 const Vector<number2> &src,
506 const number omega) const
507{
508 // check number of blocks. the sizes of the
509 // single block is checked in the function
510 // we call
511 Assert(this->n_block_cols() == 1,
512 ExcMessage("This function only works if the matrix has "
513 "a single block"));
514 Assert(this->n_block_rows() == 1,
515 ExcMessage("This function only works if the matrix has "
516 "a single block"));
517
518 // do a diagonal preconditioning. uses only
519 // the diagonal blocks of the matrix
520 this->block(0, 0).precondition_Jacobi(dst, src, omega);
521}
522
523
525
526#endif // dealii_block_sparse_matrix_h
Definition: vector.h:109
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
#define DeclException0(Exception0)
Definition: exceptions.h:464
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
#define Assert(cond, exc)
Definition: exceptions.h:1473
SmartPointer< const BlockSparsityPattern, BlockSparseMatrix< number > > sparsity_pattern
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcBlockDimensionMismatch()
static ::ExceptionBase & ExcNotQuadratic()
void print_formatted(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1.) const
static ::ExceptionBase & ExcMessage(std::string arg1)
BlockSparseMatrix(const BlockSparsityPattern &sparsity)
typename BaseClass::value_type value_type
typename BaseClass::iterator iterator
typename BaseClass::const_pointer const_pointer
BlockSparseMatrix()=default
virtual void reinit(const BlockSparsityPattern &sparsity)
size_type n_nonzero_elements() const
void precondition_Jacobi(Vector< number2 > &dst, const Vector< number2 > &src, const number omega=1.) const
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, true > > const_iterator
std::size_t memory_consumption() const
void vmult(BlockVector< block_number > &dst, const Vector< nonblock_number > &src) const
bool empty() const
void Tvmult(Vector< nonblock_number > &dst, const Vector< nonblock_number > &src) const
typename BlockType::value_type value_type
virtual ~BlockSparseMatrix() override
typename BaseClass::pointer pointer
BlockSparseMatrix & operator=(const BlockSparseMatrix &)
typename BaseClass::const_reference const_reference
size_type n_actually_nonzero_elements(const double threshold=0.0) const
void Tvmult(Vector< nonblock_number > &dst, const BlockVector< block_number > &src) const
typename BaseClass::size_type size_type
void vmult(Vector< nonblock_number > &dst, const BlockVector< block_number > &src) const
void vmult(Vector< nonblock_number > &dst, const Vector< nonblock_number > &src) const
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, false > > iterator
typename BaseClass::const_iterator const_iterator
void Tvmult(BlockVector< block_number > &dst, const BlockVector< block_number > &src) const
void Tvmult(BlockVector< block_number > &dst, const Vector< nonblock_number > &src) const
void vmult(BlockVector< block_number > &dst, const BlockVector< block_number > &src) const
typename BaseClass::BlockType BlockType
typename BaseClass::reference reference
size_type get_row_length(const size_type row) const
void precondition_Jacobi(BlockVectorType &dst, const BlockVectorType &src, const number omega=1.) const
BlockSparseMatrix & operator=(const double d)
const BlockSparsityPattern & get_sparsity_pattern() const