Reference documentation for deal.II version 9.2.0
\(\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\}}\)
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | Friends | List of all members
BlockMatrixBase< MatrixType > Class Template Reference

#include <deal.II/lac/affine_constraints.h>

Inheritance diagram for BlockMatrixBase< MatrixType >:
[legend]

Classes

struct  TemporaryData
 

Public Types

using BlockType = MatrixType
 
using value_type = typename BlockType::value_type
 
using real_type = typename numbers::NumberTraits< value_type >::real_type
 
using pointer = value_type *
 
using const_pointer = const value_type *
 
using reference = value_type &
 
using const_reference = const value_type &
 
using size_type = types::global_dof_index
 
using iterator = MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, false > >
 
using const_iterator = MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, true > >
 

Public Member Functions

 BlockMatrixBase ()=default
 
 ~BlockMatrixBase () override
 
template<class BlockMatrixType >
BlockMatrixBasecopy_from (const BlockMatrixType &source)
 
BlockTypeblock (const unsigned int row, const unsigned int column)
 
const BlockTypeblock (const unsigned int row, const unsigned int column) const
 
size_type m () const
 
size_type n () const
 
unsigned int n_block_rows () const
 
unsigned int n_block_cols () const
 
void set (const size_type i, const size_type j, const value_type value)
 
template<typename number >
void set (const std::vector< size_type > &indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=false)
 
template<typename number >
void set (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=false)
 
template<typename number >
void set (const size_type row, const std::vector< size_type > &col_indices, const std::vector< number > &values, const bool elide_zero_values=false)
 
template<typename number >
void set (const size_type row, const size_type n_cols, const size_type *col_indices, const number *values, const bool elide_zero_values=false)
 
void add (const size_type i, const size_type j, const value_type value)
 
template<typename number >
void add (const std::vector< size_type > &indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=true)
 
template<typename number >
void add (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=true)
 
template<typename number >
void add (const size_type row, const std::vector< size_type > &col_indices, const std::vector< number > &values, const bool elide_zero_values=true)
 
template<typename number >
void add (const size_type row, const size_type n_cols, const size_type *col_indices, const number *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false)
 
void add (const value_type factor, const BlockMatrixBase< MatrixType > &matrix)
 
value_type operator() (const size_type i, const size_type j) const
 
value_type el (const size_type i, const size_type j) const
 
value_type diag_element (const size_type i) const
 
void compress (::VectorOperation::values operation)
 
BlockMatrixBaseoperator*= (const value_type factor)
 
BlockMatrixBaseoperator/= (const value_type factor)
 
template<class BlockVectorType >
void vmult_add (BlockVectorType &dst, const BlockVectorType &src) const
 
template<class BlockVectorType >
void Tvmult_add (BlockVectorType &dst, const BlockVectorType &src) const
 
template<class BlockVectorType >
value_type matrix_norm_square (const BlockVectorType &v) const
 
real_type frobenius_norm () const
 
template<class BlockVectorType >
value_type matrix_scalar_product (const BlockVectorType &u, const BlockVectorType &v) const
 
template<class BlockVectorType >
value_type residual (BlockVectorType &dst, const BlockVectorType &x, const BlockVectorType &b) const
 
void print (std::ostream &out, const bool alternative_output=false) const
 
iterator begin ()
 
iterator end ()
 
iterator begin (const size_type r)
 
iterator end (const size_type r)
 
const_iterator begin () const
 
const_iterator end () const
 
const_iterator begin (const size_type r) const
 
const_iterator end (const size_type r) const
 
const BlockIndicesget_row_indices () const
 
const BlockIndicesget_column_indices () const
 
std::size_t memory_consumption () const
 

Static Public Member Functions

static ::ExceptionBaseExcIncompatibleRowNumbers (int arg1, int arg2, int arg3, int arg4)
 
static ::ExceptionBaseExcIncompatibleColNumbers (int arg1, int arg2, int arg3, int arg4)
 

Protected Member Functions

void clear ()
 
void collect_sizes ()
 
template<class BlockVectorType >
void vmult_block_block (BlockVectorType &dst, const BlockVectorType &src) const
 
template<class BlockVectorType , class VectorType >
void vmult_block_nonblock (BlockVectorType &dst, const VectorType &src) const
 
template<class BlockVectorType , class VectorType >
void vmult_nonblock_block (VectorType &dst, const BlockVectorType &src) const
 
template<class VectorType >
void vmult_nonblock_nonblock (VectorType &dst, const VectorType &src) const
 
template<class BlockVectorType >
void Tvmult_block_block (BlockVectorType &dst, const BlockVectorType &src) const
 
template<class BlockVectorType , class VectorType >
void Tvmult_block_nonblock (BlockVectorType &dst, const VectorType &src) const
 
template<class BlockVectorType , class VectorType >
void Tvmult_nonblock_block (VectorType &dst, const BlockVectorType &src) const
 
template<class VectorType >
void Tvmult_nonblock_nonblock (VectorType &dst, const VectorType &src) const
 
void prepare_add_operation ()
 
void prepare_set_operation ()
 

Protected Attributes

BlockIndices row_block_indices
 
BlockIndices column_block_indices
 
Table< 2, SmartPointer< BlockType, BlockMatrixBase< MatrixType > > > sub_objects
 

Private Attributes

TemporaryData temporary_data
 

Friends

template<typename , bool >
class BlockMatrixIterators::Accessor
 
template<typename >
class MatrixIterator
 

Detailed Description

template<typename MatrixType>
class BlockMatrixBase< MatrixType >

Blocked matrix class. The behaviour of objects of this type is almost as for the usual matrix objects, with most of the functions being implemented in both classes. The main difference is that the matrix represented by this object is composed of an array of matrices (e.g. of type SparseMatrix<number>) and all accesses to the elements of this object are relayed to accesses of the base matrices. The actual type of the individual blocks of this matrix is the type of the template argument, and can, for example be the usual SparseMatrix or PETScWrappers::SparseMatrix.

In addition to the usual matrix access and linear algebra functions, there are functions block() which allow access to the different blocks of the matrix. This may, for example, be of help when you want to implement Schur complement methods, or block preconditioners, where each block belongs to a specific component of the equation you are presently discretizing.

Note that the numbers of blocks and rows are implicitly determined by the sparsity pattern objects used.

Objects of this type are frequently used when a system of differential equations has solutions with variables that fall into different classes. For example, solutions of the Stokes or Navier-Stokes equations have dim velocity components and one pressure component. In this case, it may make sense to consider the linear system of equations as a system of 2x2 blocks, and one can construct preconditioners or solvers based on this 2x2 block structure. This class can help you in these cases, as it allows to view the matrix alternatively as one big matrix, or as a number of individual blocks.

Inheriting from this class

Since this class simply forwards its calls to the subobjects (if necessary after adjusting indices denoting which subobject is meant), this class is completely independent of the actual type of the subobject. The functions that set up block matrices and destroy them, however, have to be implemented in derived classes. These functions also have to fill the data members provided by this base class, as they are only used passively in this class.

Most of the functions take a vector or block vector argument. These functions can, in general, only successfully be compiled if the individual blocks of this matrix implement the respective functions operating on the vector type in question. For example, if you have a block sparse matrix over deal.II SparseMatrix objects, then you will likely not be able to form the matrix-vector multiplication with a block vector over PETScWrappers::SparseMatrix objects. If you attempt anyway, you will likely get a number of compiler errors.

Note
Instantiations for this template are provided for <float> and <double>; others can be generated in application programs (see the section on Template instantiations in the manual).
See also
Block (linear algebra)
Author
Wolfgang Bangerth, 2000, 2004

Definition at line 1903 of file affine_constraints.h.


The documentation for this class was generated from the following files: