Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/block_sparse_matrix.h>
Public Types | |
typedef BlockMatrixBase< SparseMatrix< number > > | BaseClass |
typedef BaseClass::BlockType | BlockType |
typedef BaseClass::value_type | value_type |
Public Types inherited from BlockMatrixBase< SparseMatrix< number > > | |
typedef SparseMatrix< number > | BlockType |
typedef BlockType::value_type | value_type |
Public Member Functions | |
Constructors and initialization | |
BlockSparseMatrix ()=default | |
BlockSparseMatrix (const BlockSparsityPattern &sparsity) | |
virtual | ~BlockSparseMatrix () |
BlockSparseMatrix & | operator= (const BlockSparseMatrix &) |
BlockSparseMatrix & | operator= (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 BlockSparsityPattern & | get_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 |
Input/Output | |
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 |
Public Member Functions inherited from BlockMatrixBase< SparseMatrix< number > > | |
BlockMatrixBase ()=default | |
~BlockMatrixBase () | |
BlockMatrixBase & | copy_from (const BlockMatrixType &source) |
BlockType & | block (const unsigned int row, const unsigned int column) |
const BlockType & | block (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) |
BlockMatrixBase & | operator*= (const value_type factor) |
BlockMatrixBase & | operator/= (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 |
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 BlockIndices & | get_row_indices () const |
const BlockIndices & | get_column_indices () const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcBlockDimensionMismatch () |
Static Public Member Functions inherited from BlockMatrixBase< SparseMatrix< number > > | |
static ::ExceptionBase & | ExcIncompatibleRowNumbers (int arg1, int arg2, int arg3, int arg4) |
static ::ExceptionBase & | ExcIncompatibleColNumbers (int arg1, int arg2, int arg3, int arg4) |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Private Attributes | |
SmartPointer< const BlockSparsityPattern, BlockSparseMatrix< number > > | sparsity_pattern |
Additional Inherited Members | |
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 |
Table< 2, SmartPointer< BlockType, BlockMatrixBase< SparseMatrix< number > > > > | sub_objects |
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.
Definition at line 50 of file block_sparse_matrix.h.
typedef BlockMatrixBase<SparseMatrix<number> > BlockSparseMatrix< number >::BaseClass |
Typedef the base class for simpler access to its own typedefs.
Definition at line 56 of file block_sparse_matrix.h.
typedef BaseClass::BlockType BlockSparseMatrix< number >::BlockType |
Typedef the type of the underlying matrix.
Definition at line 61 of file block_sparse_matrix.h.
typedef BaseClass::value_type BlockSparseMatrix< number >::value_type |
Import the typedefs from the base class.
Definition at line 66 of file block_sparse_matrix.h.
|
default |
Constructor; initializes the matrix to be empty, without any structure, i.e. the matrix is not usable at all. This constructor is therefore only useful for matrices which are members of a class. All other matrices should be created at a point in the data flow where all necessary information is available.
You have to initialize the matrix before usage with reinit(BlockSparsityPattern). The number of blocks per row and column are then determined by that function.
BlockSparseMatrix< number >::BlockSparseMatrix | ( | const BlockSparsityPattern & | sparsity | ) |
Constructor. Takes the given matrix sparsity structure to represent the sparsity pattern of this matrix. You can change the sparsity pattern later on by calling the reinit() function.
This constructor initializes all sub-matrices with the sub-sparsity pattern within the argument.
You have to make sure that the lifetime of the sparsity structure is at least as long as that of this matrix or as long as reinit() is not called with a new sparsity structure.
|
virtual |
Destructor.
BlockSparseMatrix& BlockSparseMatrix< number >::operator= | ( | const BlockSparseMatrix< number > & | ) |
Pseudo copy operator only copying empty objects. The sizes of the block matrices need to be the same.
|
inline |
This operator assigns a scalar to a matrix. Since this does usually not make much sense (should we set all matrix entries to this value? Only the nonzero entries of the sparsity pattern?), this operation is only allowed if the actual value to be assigned is zero. This operator only exists to allow for the obvious notation matrix=0
, which sets all elements of the matrix to zero, but keep the sparsity pattern previously used.
Definition at line 367 of file block_sparse_matrix.h.
void BlockSparseMatrix< number >::clear | ( | ) |
Release all memory and return to a state just like after having called the default constructor. It also forgets the sparsity pattern it was previously tied to.
This calls SparseMatrix::clear on all sub-matrices and then resets this object to have no blocks at all.
|
virtual |
Reinitialize the sparse matrix with the given sparsity pattern. The latter tells the matrix how many nonzero elements there need to be reserved.
Basically, this function only calls SparseMatrix::reinit() of the sub- matrices with the block sparsity patterns of the parameter.
You have to make sure that the lifetime of the sparsity structure is at least as long as that of this matrix or as long as reinit(const SparsityPattern &) is not called with a new sparsity structure.
The elements of the matrix are set to zero by this function.
bool BlockSparseMatrix< number >::empty | ( | ) | const |
Return whether the object is empty. It is empty if either both dimensions are zero or no BlockSparsityPattern is associated.
size_type BlockSparseMatrix< number >::get_row_length | ( | const size_type | row | ) | const |
Return the number of entries in a specific row.
size_type BlockSparseMatrix< number >::n_nonzero_elements | ( | ) | const |
Return the number of nonzero elements of this matrix. Actually, it returns the number of entries in the sparsity pattern; if any of the entries should happen to be zero, it is counted anyway.
size_type BlockSparseMatrix< number >::n_actually_nonzero_elements | ( | const double | threshold = 0.0 | ) | const |
Return the number of actually nonzero elements. Just counts the number of actually nonzero elements (with absolute value larger than threshold) of all the blocks.
const BlockSparsityPattern& BlockSparseMatrix< number >::get_sparsity_pattern | ( | ) | const |
Return a (constant) reference to the underlying sparsity pattern of this matrix.
Though the return value is declared const
, you should be aware that it may change if you call any nonconstant function of objects which operate on it.
std::size_t BlockSparseMatrix< number >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
|
inline |
Matrix-vector multiplication: let \(dst = M*src\) with \(M\) being this matrix.
Definition at line 384 of file block_sparse_matrix.h.
|
inline |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block column.
Definition at line 397 of file block_sparse_matrix.h.
|
inline |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block row.
Definition at line 410 of file block_sparse_matrix.h.
|
inline |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block.
Definition at line 422 of file block_sparse_matrix.h.
|
inline |
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.
Definition at line 434 of file block_sparse_matrix.h.
|
inline |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block row.
Definition at line 447 of file block_sparse_matrix.h.
|
inline |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block column.
Definition at line 460 of file block_sparse_matrix.h.
|
inline |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block.
Definition at line 472 of file block_sparse_matrix.h.
|
inline |
Apply the Jacobi preconditioner, which multiplies every element of the src
vector by the inverse of the respective diagonal element and multiplies the result with the relaxation parameter omega
.
All diagonal blocks must be square matrices for this operation.
Definition at line 485 of file block_sparse_matrix.h.
|
inline |
Apply the Jacobi preconditioner to a simple vector.
The matrix must be a single square block for this.
Definition at line 510 of file block_sparse_matrix.h.
void BlockSparseMatrix< number >::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 |
Print the matrix in the usual format, i.e. as a matrix and not as a list of nonzero elements. For better readability, elements not in the matrix are displayed as empty space, while matrix elements which are explicitly set to zero are displayed as such.
The parameters allow for a flexible setting of the output format: precision
and scientific
are used to determine the number format, where scientific = false
means fixed point notation. A zero entry for width
makes the function compute a width, but it may be changed to a positive value, if output is crude.
Additionally, a character for an empty value may be specified.
Finally, the whole matrix can be multiplied with a common denominator to produce more readable output, even integers.
|
private |
Pointer to the block sparsity pattern used for this matrix. In order to guarantee that it is not deleted while still in use, we subscribe to it using the SmartPointer class.
Definition at line 354 of file block_sparse_matrix.h.