Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/matrix_block.h>
Public Types | |
using | size_type = types::global_dof_index |
using | value_type = typename MatrixType::value_type |
Public Member Functions | |
MatrixBlock () | |
MatrixBlock (const MatrixBlock< MatrixType > &M) | |
MatrixBlock (size_type i, size_type j) | |
void | reinit (const BlockSparsityPattern &sparsity) |
void | add (const size_type i, const size_type j, const typename MatrixType::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_index, 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) |
template<class VectorType > | |
void | vmult (VectorType &w, const VectorType &v) const |
template<class VectorType > | |
void | vmult_add (VectorType &w, const VectorType &v) const |
template<class VectorType > | |
void | Tvmult (VectorType &w, const VectorType &v) const |
template<class VectorType > | |
void | Tvmult_add (VectorType &w, const VectorType &v) 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 (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 Public Member Functions | |
static ::ExceptionBase & | ExcBlockIndexMismatch (size_type arg1, size_type arg2) |
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) |
Public Attributes | |
size_type | row |
size_type | column |
MatrixType | matrix |
Private Attributes | |
BlockIndices | row_indices |
BlockIndices | column_indices |
A wrapper around a matrix object, storing the coordinates in a block matrix as well.
This class is an alternative to BlockMatrixBase, if you only want to generate a single block of the system, not the whole system. Using the add() functions of this class, it is possible to use the standard assembling functions used for block matrices, but only enter in one of the blocks and still avoiding the index computations involved. The reason for this class is, that we may need a different number of matrices for different blocks in a block system. For example, a preconditioner for the Oseen system can be built as a block system, where the pressure block is of the form M-1FA-1 with M the pressure mass matrix, A the pressure Laplacian and F the advection diffusion operator applied to the pressure space. Since only a single matrix is needed for the other blocks, using BlockSparseMatrix or similar would be a waste of memory.
While the add() functions make a MatrixBlock appear like a block matrix for assembling, the functions vmult(), Tvmult(), vmult_add(), and Tvmult_add() make it behave like a MatrixType, when it comes to applying it to a vector. This behavior allows us to store MatrixBlock objects in vectors, for instance in MGLevelObject without extracting the matrix first.
MatrixBlock comes handy when using BlockMatrixArray. Once the MatrixBlock has been properly initialized and filled, it can be used in the simplest case as:
Here, we have not gained very much, except that we do not need to set up empty blocks in the block system.
Definition at line 38 of file matrix_block.h.
using MatrixBlock< MatrixType >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 115 of file matrix_block.h.
using MatrixBlock< MatrixType >::value_type = typename MatrixType::value_type |
Declare a type for matrix entries.
Definition at line 120 of file matrix_block.h.
|
inline |
Constructor rendering an uninitialized object.
Definition at line 632 of file matrix_block.h.
|
inline |
Copy constructor.
Definition at line 639 of file matrix_block.h.
|
inline |
Constructor setting block coordinates, but not initializing the matrix.
Definition at line 650 of file matrix_block.h.
|
inline |
Reinitialize the matrix for a new BlockSparsityPattern. This adjusts the matrix as well as the row_indices and column_indices.
Definition at line 658 of file matrix_block.h.
|
inline |
Add value
to the element (i,j). Throws an error if the entry does not exist or if it is in a different block.
Definition at line 680 of file matrix_block.h.
|
inline |
Add all elements in a FullMatrix into sparse matrix locations given by indices
. This function assumes a quadratic sparse matrix and a quadratic full_matrix. The global locations are translated into locations in this block and ExcBlockIndexMismatch is thrown, if the global index does not point into the block referred to by row and column.
elide_zero_values
is currently ignored.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.
Definition at line 763 of file matrix_block.h.
|
inline |
Add all elements in a FullMatrix into global locations given by row_indices
and col_indices
, respectively. The global locations are translated into locations in this block and ExcBlockIndexMismatch is thrown, if the global index does not point into the block referred to by row and column.
elide_zero_values
is currently ignored.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.
Definition at line 701 of file matrix_block.h.
|
inline |
Set several elements in the specified row of the matrix with column indices as given by col_indices
to the respective value. This is the function doing the actual work for the ones adding full matrices. The global locations row_index
and col_indices
are translated into locations in this block and ExcBlockIndexMismatch is thrown, if the global index does not point into the block referred to by row and column.
elide_zero_values
is currently ignored.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.
Definition at line 786 of file matrix_block.h.
|
inline |
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.
Definition at line 724 of file matrix_block.h.
|
inline |
Matrix-vector-multiplication, forwarding to the same function in MatrixType. No index computations are done, thus, the vectors need to have sizes matching matrix.
Definition at line 806 of file matrix_block.h.
|
inline |
Matrix-vector-multiplication, forwarding to the same function in MatrixType. No index computations are done, thus, the vectors need to have sizes matching matrix.
Definition at line 815 of file matrix_block.h.
|
inline |
Matrix-vector-multiplication, forwarding to the same function in MatrixType. No index computations are done, thus, the vectors need to have sizes matching matrix.
Definition at line 824 of file matrix_block.h.
|
inline |
Matrix-vector-multiplication, forwarding to the same function in MatrixType. No index computations are done, thus, the vectors need to have sizes matching matrix.
Definition at line 833 of file matrix_block.h.
|
inline |
The memory used by this object.
Definition at line 841 of file matrix_block.h.
size_type MatrixBlock< MatrixType >::row |
Row coordinate. This is the position of the data member matrix on the global matrix.
Definition at line 298 of file matrix_block.h.
size_type MatrixBlock< MatrixType >::column |
Column coordinate. This is the position of the data member matrix on the global matrix.
Definition at line 303 of file matrix_block.h.
MatrixType MatrixBlock< MatrixType >::matrix |
The matrix itself
Definition at line 308 of file matrix_block.h.
|
private |
The row BlockIndices of the whole system. Using row(), this allows us to find the index of the first row degree of freedom for this block.
Definition at line 315 of file matrix_block.h.
|
private |
The column BlockIndices of the whole system. Using column(), this allows us to find the index of the first column degree of freedom for this block.
Definition at line 321 of file matrix_block.h.