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
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | 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
 
using real_type = typename numbers::NumberTraits< value_type >::real_type
 

Public Member Functions

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
 
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
 

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 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

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

Private Attributes

TemporaryData temporary_data
 

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 ()
 

Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 
static std::mutex mutex
 
void check_no_subscribers () const noexcept
 

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)

Definition at line 49 of file block_sparse_matrix.h.

Member Typedef Documentation

◆ real_type

using BlockMatrixBase< SparseMatrix< number > >::real_type = typename numbers::NumberTraits<value_type>::real_type
inherited

Definition at line 362 of file block_matrix_base.h.

Member Function Documentation

◆ copy_from()

BlockMatrixBase & BlockMatrixBase< SparseMatrix< number > >::copy_from ( const BlockMatrixType &  source)
inherited

Copy the matrix given as argument into the current object.

Copying matrices is an expensive operation that we do not want to happen by accident through compiler generated code for operator=. (This would happen, for example, if one accidentally declared a function argument of the current type by value rather than by reference.) The functionality of copying matrices is implemented in this member function instead. All copy operations of objects of this type therefore require an explicit function call.

The source matrix may be a matrix of arbitrary type, as long as its data type is convertible to the data type of this matrix.

The function returns a reference to this.

◆ block() [1/2]

BlockType & BlockMatrixBase< SparseMatrix< number > >::block ( const unsigned int  row,
const unsigned int  column 
)
inherited

Access the block with the given coordinates.

◆ block() [2/2]

const BlockType & BlockMatrixBase< SparseMatrix< number > >::block ( const unsigned int  row,
const unsigned int  column 
) const
inherited

Access the block with the given coordinates. Version for constant objects.

◆ m()

size_type BlockMatrixBase< SparseMatrix< number > >::m ( ) const
inherited

Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).

◆ n()

size_type BlockMatrixBase< SparseMatrix< number > >::n ( ) const
inherited

Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).

◆ n_block_rows()

unsigned int BlockMatrixBase< SparseMatrix< number > >::n_block_rows ( ) const
inherited

Return the number of blocks in a column. Returns zero if no sparsity pattern is presently associated to this matrix.

◆ n_block_cols()

unsigned int BlockMatrixBase< SparseMatrix< number > >::n_block_cols ( ) const
inherited

Return the number of blocks in a row. Returns zero if no sparsity pattern is presently associated to this matrix.

◆ set() [1/5]

void BlockMatrixBase< SparseMatrix< number > >::set ( const size_type  i,
const size_type  j,
const value_type  value 
)
inherited

Set the element (i,j) to value. Throws an error if the entry does not exist or if value is not a finite number. Still, it is allowed to store zero values in non-existent fields.

◆ set() [2/5]

void BlockMatrixBase< SparseMatrix< number > >::set ( const std::vector< size_type > &  indices,
const FullMatrix< number > &  full_matrix,
const bool  elide_zero_values = false 
)
inherited

Set all elements given in a FullMatrix into the sparse matrix locations given by indices. In other words, this function writes the elements in full_matrix into the calling matrix, using the local-to-global indexing specified by indices for both the rows and the columns of the matrix. This function assumes a quadratic sparse matrix and a quadratic full_matrix, the usual situation in FE calculations.

The optional parameter elide_zero_values can be used to specify whether zero values should be set anyway or they should be filtered away (and not change the previous content in the respective element if it exists). The default value is false, i.e., even zero values are treated.

◆ set() [3/5]

void BlockMatrixBase< SparseMatrix< number > >::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 
)
inherited

Same function as before, but now including the possibility to use rectangular full_matrices and different local-to-global indexing on rows and columns, respectively.

◆ set() [4/5]

void BlockMatrixBase< SparseMatrix< number > >::set ( const size_type  row,
const std::vector< size_type > &  col_indices,
const std::vector< number > &  values,
const bool  elide_zero_values = false 
)
inherited

Set several elements in the specified row of the matrix with column indices as given by col_indices to the respective value.

The optional parameter elide_zero_values can be used to specify whether zero values should be set anyway or they should be filtered away (and not change the previous content in the respective element if it exists). The default value is false, i.e., even zero values are treated.

◆ set() [5/5]

void BlockMatrixBase< SparseMatrix< number > >::set ( const size_type  row,
const size_type  n_cols,
const size_type col_indices,
const number *  values,
const bool  elide_zero_values = false 
)
inherited

Set several elements to values given by values in a given row in columns given by col_indices into the sparse matrix.

The optional parameter elide_zero_values can be used to specify whether zero values should be inserted anyway or they should be filtered away. The default value is false, i.e., even zero values are inserted/replaced.

◆ add() [1/6]

void BlockMatrixBase< SparseMatrix< number > >::add ( const size_type  i,
const size_type  j,
const value_type  value 
)
inherited

Add value to the element (i,j). Throws an error if the entry does not exist or if value is not a finite number. Still, it is allowed to store zero values in non-existent fields.

◆ add() [2/6]

void BlockMatrixBase< SparseMatrix< number > >::add ( const std::vector< size_type > &  indices,
const FullMatrix< number > &  full_matrix,
const bool  elide_zero_values = true 
)
inherited

Add all elements given in a FullMatrix<double> into sparse matrix locations given by indices. In other words, this function adds the elements in full_matrix to the respective entries in calling matrix, using the local-to-global indexing specified by indices for both the rows and the columns of the matrix. This function assumes a quadratic sparse matrix and a quadratic full_matrix, the usual situation in FE calculations.

The optional parameter elide_zero_values can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true, i.e., zero values won't be added into the matrix.

◆ add() [3/6]

void BlockMatrixBase< SparseMatrix< number > >::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 
)
inherited

Same function as before, but now including the possibility to use rectangular full_matrices and different local-to-global indexing on rows and columns, respectively.

◆ add() [4/6]

void BlockMatrixBase< SparseMatrix< number > >::add ( const size_type  row,
const std::vector< size_type > &  col_indices,
const std::vector< number > &  values,
const bool  elide_zero_values = true 
)
inherited

Set several elements in the specified row of the matrix with column indices as given by col_indices to the respective value.

The optional parameter elide_zero_values can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true, i.e., zero values won't be added into the matrix.

◆ add() [5/6]

void BlockMatrixBase< SparseMatrix< number > >::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 
)
inherited

Add an array of values given by values in the given global matrix row at columns specified by col_indices in the sparse matrix.

The optional parameter elide_zero_values can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true, i.e., zero values won't be added into the matrix.

◆ add() [6/6]

void BlockMatrixBase< SparseMatrix< number > >::add ( const value_type  factor,
const BlockMatrixBase< SparseMatrix< number > > &  matrix 
)
inherited

Add matrix scaled by factor to this matrix, i.e. the matrix factor*matrix is added to this. If the sparsity pattern of the calling matrix does not contain all the elements in the sparsity pattern of the input matrix, this function will throw an exception.

Depending on MatrixType, however, additional restrictions might arise. Some sparse matrix formats require matrix to be based on the same sparsity pattern as the calling matrix.

◆ operator()()

value_type BlockMatrixBase< SparseMatrix< number > >::operator() ( const size_type  i,
const size_type  j 
) const
inherited

Return the value of the entry (i,j). This may be an expensive operation and you should always take care where to call this function. In order to avoid abuse, this function throws an exception if the wanted element does not exist in the matrix.

◆ el()

value_type BlockMatrixBase< SparseMatrix< number > >::el ( const size_type  i,
const size_type  j 
) const
inherited

This function is mostly like operator()() in that it returns the value of the matrix entry (i,j). The only difference is that if this entry does not exist in the sparsity pattern, then instead of raising an exception, zero is returned. While this may be convenient in some cases, note that it is simple to write algorithms that are slow compared to an optimal solution, since the sparsity of the matrix is not used.

◆ diag_element()

value_type BlockMatrixBase< SparseMatrix< number > >::diag_element ( const size_type  i) const
inherited

Return the main diagonal element in the ith row. This function throws an error if the matrix is not quadratic and also if the diagonal blocks of the matrix are not quadratic.

This function is considerably faster than the operator()(), since for quadratic matrices, the diagonal entry may be the first to be stored in each row and access therefore does not involve searching for the right column number.

◆ compress()

void BlockMatrixBase< SparseMatrix< number > >::compress ( ::VectorOperation::values  operation)
inherited

Call the compress() function on all the subblocks of the matrix.

See Compressing distributed objects for more information.

◆ operator*=()

BlockMatrixBase & BlockMatrixBase< SparseMatrix< number > >::operator*= ( const value_type  factor)
inherited

Multiply the entire matrix by a fixed factor.

◆ operator/=()

BlockMatrixBase & BlockMatrixBase< SparseMatrix< number > >::operator/= ( const value_type  factor)
inherited

Divide the entire matrix by a fixed factor.

◆ vmult_add()

void BlockMatrixBase< SparseMatrix< number > >::vmult_add ( BlockVectorType &  dst,
const BlockVectorType &  src 
) const
inherited

Adding Matrix-vector multiplication. Add \(M*src\) on \(dst\) with \(M\) being this matrix.

◆ Tvmult_add()

void BlockMatrixBase< SparseMatrix< number > >::Tvmult_add ( BlockVectorType &  dst,
const BlockVectorType &  src 
) const
inherited

Adding Matrix-vector multiplication. Add MTsrc to dst with M being this matrix. This function does the same as vmult_add() but takes the transposed matrix.

◆ matrix_norm_square()

value_type BlockMatrixBase< SparseMatrix< number > >::matrix_norm_square ( const BlockVectorType &  v) const
inherited

Return the norm of the vector v with respect to the norm induced by this matrix, i.e. vTMv). This is useful, e.g. in the finite element context, where the LT-norm of a function equals the matrix norm with respect to the mass matrix of the vector representing the nodal values of the finite element function. Note that even though the function's name might suggest something different, for historic reasons not the norm but its square is returned, as defined above by the scalar product.

Obviously, the matrix needs to be square for this operation.

◆ frobenius_norm()

real_type BlockMatrixBase< SparseMatrix< number > >::frobenius_norm ( ) const
inherited

Return the frobenius norm of the matrix, i.e. the square root of the sum of squares of all entries in the matrix.

◆ matrix_scalar_product()

value_type BlockMatrixBase< SparseMatrix< number > >::matrix_scalar_product ( const BlockVectorType &  u,
const BlockVectorType &  v 
) const
inherited

Compute the matrix scalar product \(\left(u,Mv\right)\).

◆ residual()

value_type BlockMatrixBase< SparseMatrix< number > >::residual ( BlockVectorType &  dst,
const BlockVectorType &  x,
const BlockVectorType &  b 
) const
inherited

Compute the residual r=b-Ax. Write the residual into dst.

◆ print()

void BlockMatrixBase< SparseMatrix< number > >::print ( std::ostream &  out,
const bool  alternative_output = false 
) const
inherited

Print the matrix to the given stream, using the format (line,col) value, i.e. one nonzero entry of the matrix per line. The optional flag outputs the sparsity pattern in a different style according to the underlying sparse matrix type.

◆ begin() [1/4]

iterator BlockMatrixBase< SparseMatrix< number > >::begin ( )
inherited

Iterator starting at the first entry.

◆ begin() [2/4]

iterator BlockMatrixBase< SparseMatrix< number > >::begin ( const size_type  r)
inherited

Iterator starting at the first entry of row r.

◆ begin() [3/4]

const_iterator BlockMatrixBase< SparseMatrix< number > >::begin ( ) const
inherited

Iterator starting at the first entry.

◆ begin() [4/4]

const_iterator BlockMatrixBase< SparseMatrix< number > >::begin ( const size_type  r) const
inherited

Iterator starting at the first entry of row r.

◆ end() [1/4]

iterator BlockMatrixBase< SparseMatrix< number > >::end ( )
inherited

Final iterator.

◆ end() [2/4]

iterator BlockMatrixBase< SparseMatrix< number > >::end ( const size_type  r)
inherited

Final iterator of row r.

◆ end() [3/4]

const_iterator BlockMatrixBase< SparseMatrix< number > >::end ( ) const
inherited

Final iterator.

◆ end() [4/4]

const_iterator BlockMatrixBase< SparseMatrix< number > >::end ( const size_type  r) const
inherited

Final iterator of row r.

◆ get_row_indices()

const BlockIndices & BlockMatrixBase< SparseMatrix< number > >::get_row_indices ( ) const
inherited

Return a reference to the underlying BlockIndices data of the rows.

◆ get_column_indices()

const BlockIndices & BlockMatrixBase< SparseMatrix< number > >::get_column_indices ( ) const
inherited

Return a reference to the underlying BlockIndices data of the columns.

◆ ExcIncompatibleRowNumbers()

static ::ExceptionBase & BlockMatrixBase< SparseMatrix< number > >::ExcIncompatibleRowNumbers ( int  arg1,
int  arg2,
int  arg3,
int  arg4 
)
staticinherited

Exception

Note
The message that will be printed by this exception reads:
<< "The blocks [" << arg1 << ',' << arg2 << "] and [" << arg3 << ',' << arg4 << "] have differing row numbers."

◆ ExcIncompatibleColNumbers()

static ::ExceptionBase & BlockMatrixBase< SparseMatrix< number > >::ExcIncompatibleColNumbers ( int  arg1,
int  arg2,
int  arg3,
int  arg4 
)
staticinherited

Exception

Note
The message that will be printed by this exception reads:
<< "The blocks [" << arg1 << ',' << arg2 << "] and [" << arg3 << ',' << arg4 << "] have differing column numbers."

◆ collect_sizes()

void BlockMatrixBase< SparseMatrix< number > >::collect_sizes ( )
protectedinherited

This function collects the sizes of the sub-objects and stores them in internal arrays, in order to be able to relay global indices into the matrix to indices into the subobjects. You must call this function each time after you have changed the size of the sub-objects.

Derived classes should call this function whenever the size of the sub- objects has changed and the X_block_indices arrays need to be updated.

Note that this function is not public since not all derived classes need to export its interface. For example, for the usual deal.II SparseMatrix class, the sizes are implicitly determined whenever reinit() is called, and individual blocks cannot be resized. For that class, this function therefore does not have to be public. On the other hand, for the PETSc classes, there is no associated sparsity pattern object that determines the block sizes, and for these the function needs to be publicly available. These classes therefore export this function.

◆ vmult_block_block()

void BlockMatrixBase< SparseMatrix< number > >::vmult_block_block ( BlockVectorType &  dst,
const BlockVectorType &  src 
) const
protectedinherited

Matrix-vector multiplication: let \(dst = M*src\) with \(M\) being this matrix.

Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.

◆ vmult_block_nonblock()

void BlockMatrixBase< SparseMatrix< number > >::vmult_block_nonblock ( BlockVectorType &  dst,
const VectorType &  src 
) const
protectedinherited

Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block column.

Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.

◆ vmult_nonblock_block()

void BlockMatrixBase< SparseMatrix< number > >::vmult_nonblock_block ( VectorType &  dst,
const BlockVectorType &  src 
) const
protectedinherited

Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block row.

Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.

◆ vmult_nonblock_nonblock()

void BlockMatrixBase< SparseMatrix< number > >::vmult_nonblock_nonblock ( VectorType &  dst,
const VectorType &  src 
) const
protectedinherited

Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block.

Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.

◆ Tvmult_block_block()

void BlockMatrixBase< SparseMatrix< number > >::Tvmult_block_block ( BlockVectorType &  dst,
const BlockVectorType &  src 
) const
protectedinherited

Matrix-vector multiplication: let \(dst = M^T*src\) with \(M\) being this matrix. This function does the same as vmult() but takes the transposed matrix.

Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.

◆ Tvmult_block_nonblock()

void BlockMatrixBase< SparseMatrix< number > >::Tvmult_block_nonblock ( BlockVectorType &  dst,
const VectorType &  src 
) const
protectedinherited

Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block row.

Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.

◆ Tvmult_nonblock_block()

void BlockMatrixBase< SparseMatrix< number > >::Tvmult_nonblock_block ( VectorType &  dst,
const BlockVectorType &  src 
) const
protectedinherited

Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block column.

Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.

◆ Tvmult_nonblock_nonblock()

void BlockMatrixBase< SparseMatrix< number > >::Tvmult_nonblock_nonblock ( VectorType &  dst,
const VectorType &  src 
) const
protectedinherited

Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block.

Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.

◆ prepare_add_operation()

void BlockMatrixBase< SparseMatrix< number > >::prepare_add_operation ( )
protectedinherited

Some matrix types, in particular PETSc, need to synchronize set and add operations. This has to be done for all matrices in the BlockMatrix. This routine prepares adding of elements by notifying all blocks. Called by all internal routines before adding elements.

◆ prepare_set_operation()

void BlockMatrixBase< SparseMatrix< number > >::prepare_set_operation ( )
protectedinherited

Notifies all blocks to let them prepare for setting elements, see prepare_add_operation().

Member Data Documentation

◆ row_block_indices

BlockIndices BlockMatrixBase< SparseMatrix< number > >::row_block_indices
protectedinherited

Index arrays for rows and columns.

Definition at line 846 of file block_matrix_base.h.

◆ column_block_indices

BlockIndices BlockMatrixBase< SparseMatrix< number > >::column_block_indices
protectedinherited

Definition at line 847 of file block_matrix_base.h.

◆ sub_objects

Table<2, SmartPointer<BlockType, BlockMatrixBase<SparseMatrix< number > > > > BlockMatrixBase< SparseMatrix< number > >::sub_objects
protectedinherited

Array of sub-matrices.

Definition at line 852 of file block_matrix_base.h.

◆ temporary_data

TemporaryData BlockMatrixBase< SparseMatrix< number > >::temporary_data
privateinherited

A set of scratch arrays that can be used by the add() and set() functions that take pointers to data to pre-sort indices before use. Access from multiple threads is synchronized via the mutex variable that is part of the structure.

Definition at line 1064 of file block_matrix_base.h.


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