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\}}\)
Public Types | List of all members
BlockSparseMatrix< number > Class Template Reference

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

Inheritance diagram for BlockSparseMatrix< number >:
[legend]

Public Types

using BaseClass = BlockMatrixBase< SparseMatrix< number > >
 
using BlockType = typename BaseClass::BlockType
 
using value_type = typename BaseClass::value_type
 
using pointer = typename BaseClass::pointer
 
using const_pointer = typename BaseClass::const_pointer
 
using reference = typename BaseClass::reference
 
using const_reference = typename BaseClass::const_reference
 
using size_type = typename BaseClass::size_type
 
using iterator = typename BaseClass::iterator
 
using const_iterator = typename BaseClass::const_iterator
 
- Public Types inherited from BlockMatrixBase< SparseMatrix< number > >
using BlockType = SparseMatrix< number >
 
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

Constructors and initialization
 BlockSparseMatrix ()=default
 
 BlockSparseMatrix (const BlockSparsityPattern &sparsity)
 
virtual ~BlockSparseMatrix () override
 
BlockSparseMatrixoperator= (const BlockSparseMatrix &)
 
BlockSparseMatrixoperator= (const double d)
 
void clear ()
 
virtual void reinit (const BlockSparsityPattern &sparsity)
 
Information on the matrix
bool empty () const
 
size_type get_row_length (const size_type row) const
 
size_type n_nonzero_elements () const
 
size_type n_actually_nonzero_elements (const double threshold=0.0) const
 
const BlockSparsityPatternget_sparsity_pattern () const
 
std::size_t memory_consumption () const
 
Multiplications
template<typename block_number >
void vmult (BlockVector< block_number > &dst, const BlockVector< block_number > &src) const
 
template<typename block_number , typename nonblock_number >
void vmult (BlockVector< block_number > &dst, const Vector< nonblock_number > &src) const
 
template<typename block_number , typename nonblock_number >
void vmult (Vector< nonblock_number > &dst, const BlockVector< block_number > &src) const
 
template<typename nonblock_number >
void vmult (Vector< nonblock_number > &dst, const Vector< nonblock_number > &src) const
 
template<typename block_number >
void Tvmult (BlockVector< block_number > &dst, const BlockVector< block_number > &src) const
 
template<typename block_number , typename nonblock_number >
void Tvmult (BlockVector< block_number > &dst, const Vector< nonblock_number > &src) const
 
template<typename block_number , typename nonblock_number >
void Tvmult (Vector< nonblock_number > &dst, const BlockVector< block_number > &src) const
 
template<typename nonblock_number >
void Tvmult (Vector< nonblock_number > &dst, const Vector< nonblock_number > &src) const
 
Preconditioning methods
template<class BlockVectorType >
void precondition_Jacobi (BlockVectorType &dst, const BlockVectorType &src, const number omega=1.) const
 
template<typename number2 >
void precondition_Jacobi (Vector< number2 > &dst, const Vector< number2 > &src, const number omega=1.) const
 
- Public Member Functions inherited from BlockMatrixBase< SparseMatrix< number > >
 BlockMatrixBase ()=default
 
 ~BlockMatrixBase () override
 
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)
 
void set (const std::vector< size_type > &indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=false)
 
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)
 
void set (const size_type row, const std::vector< size_type > &col_indices, const std::vector< number > &values, const bool elide_zero_values=false)
 
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)
 
void add (const std::vector< size_type > &indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=true)
 
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)
 
void add (const size_type row, const std::vector< size_type > &col_indices, const std::vector< number > &values, const bool elide_zero_values=true)
 
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< SparseMatrix< number > > &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)
 
void vmult_add (BlockVectorType &dst, const BlockVectorType &src) const
 
void Tvmult_add (BlockVectorType &dst, const BlockVectorType &src) const
 
value_type matrix_norm_square (const BlockVectorType &v) const
 
real_type frobenius_norm () const
 
value_type matrix_scalar_product (const BlockVectorType &u, const BlockVectorType &v) const
 
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 begin (const size_type r)
 
const_iterator begin () const
 
const_iterator begin (const size_type r) const
 
iterator end ()
 
iterator end (const size_type r)
 
const_iterator end () 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
 

Input/Output

SmartPointer< const BlockSparsityPattern, BlockSparseMatrix< number > > sparsity_pattern
 
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 ::ExceptionBaseExcBlockDimensionMismatch ()
 

Additional Inherited Members

- Static Public Member Functions inherited from BlockMatrixBase< SparseMatrix< number > >
static ::ExceptionBaseExcIncompatibleRowNumbers (int arg1, int arg2, int arg3, int arg4)
 
static ::ExceptionBaseExcIncompatibleColNumbers (int arg1, int arg2, int arg3, int arg4)
 
- Protected Member Functions inherited from BlockMatrixBase< SparseMatrix< number > >
void clear ()
 
void collect_sizes ()
 
void vmult_block_block (BlockVectorType &dst, const BlockVectorType &src) const
 
void vmult_block_nonblock (BlockVectorType &dst, const VectorType &src) const
 
void vmult_nonblock_block (VectorType &dst, const BlockVectorType &src) const
 
void vmult_nonblock_nonblock (VectorType &dst, const VectorType &src) const
 
void Tvmult_block_block (BlockVectorType &dst, const BlockVectorType &src) const
 
void Tvmult_block_nonblock (BlockVectorType &dst, const VectorType &src) const
 
void Tvmult_nonblock_block (VectorType &dst, const BlockVectorType &src) const
 
void Tvmult_nonblock_nonblock (VectorType &dst, const VectorType &src) const
 
void prepare_add_operation ()
 
void prepare_set_operation ()
 
- Protected Attributes inherited from BlockMatrixBase< SparseMatrix< number > >
BlockIndices row_block_indices
 
BlockIndices column_block_indices
 
Table< 2, SmartPointer< BlockType, BlockMatrixBase< SparseMatrix< number > > > > sub_objects
 

Detailed Description

template<typename number>
class BlockSparseMatrix< number >

Blocked sparse matrix based on the SparseMatrix class. This class implements the functions that are specific to the SparseMatrix base objects for a blocked sparse matrix, and leaves the actual work relaying most of the calls to the individual blocks to the functions implemented in the base class. See there also for a description of when this class is useful.

See also
Block (linear algebra)
Author
Wolfgang Bangerth, 2000, 2004

Definition at line 50 of file block_sparse_matrix.h.


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