Matrix and Vector classes that support a block recursive structure capable of representing the natural structure from Finite Element discretisations.
More...
|
file | matrixmatrix.hh |
| provides functions for sparse matrix matrix multiplication.
|
|
file | matrixutils.hh |
| Some handy generic functions for ISTL matrices.
|
|
|
struct | Dune::CompressionStatistics< size_type > |
| Statistics about compression achieved in implicit mode. More...
|
|
class | Dune::ImplicitMatrixBuilder< M_ > |
| A wrapper for uniform access to the BCRSMatrix during and after the build stage in implicit build mode. More...
|
|
class | Dune::BCRSMatrix< B, A > |
| A sparse block matrix with compressed row storage. More...
|
|
class | Dune::BDMatrix< B, A > |
| A block-diagonal matrix. More...
|
|
class | Dune::BTDMatrix< B, A > |
| A block-tridiagonal matrix. More...
|
|
class | Dune::BlockVector< B, A > |
| A vector of blocks with memory management. More...
|
|
class | Dune::Matrix< T, A > |
| A generic dynamic dense matrix. More...
|
|
struct | Dune::MatMultMatResult< M1, M2 > |
| Helper TMP to get the result type of a sparse matrix matrix multiplication ( \(C=A*B\)) More...
|
|
struct | Dune::TransposedMatMultMatResult< M1, M2 > |
| Helper TMP to get the result type of a sparse matrix matrix multiplication ( \(C=A*B\)) More...
|
|
struct | Dune::CheckIfDiagonalPresent< Matrix, blocklevel, l > |
| Check whether the a matrix has diagonal values on blocklevel recursion levels. More...
|
|
class | Dune::VariableBlockVector< B, A > |
| A Vector of blocks with different blocksizes. More...
|
|
|
template<class T , class A , class A1 , class A2 , int n, int m, int k> |
void | Dune::matMultTransposeMat (BCRSMatrix< FieldMatrix< T, n, k >, A > &res, const BCRSMatrix< FieldMatrix< T, n, m >, A1 > &mat, const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &matt, bool tryHard=false) |
| Calculate product of a sparse matrix with a transposed sparse matrices ( \(C=A*B^T\)). More...
|
|
template<class T , class A , class A1 , class A2 , int n, int m, int k> |
void | Dune::matMultMat (BCRSMatrix< FieldMatrix< T, n, m >, A > &res, const BCRSMatrix< FieldMatrix< T, n, k >, A1 > &mat, const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &matt, bool tryHard=false) |
| Calculate product of two sparse matrices ( \(C=A*B\)). More...
|
|
template<class T , class A , class A1 , class A2 , int n, int m, int k> |
void | Dune::transposeMatMultMat (BCRSMatrix< FieldMatrix< T, n, m >, A > &res, const BCRSMatrix< FieldMatrix< T, k, n >, A1 > &mat, const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &matt, bool tryHard=false) |
| Calculate product of a transposed sparse matrix with another sparse matrices ( \(C=A^T*B\)). More...
|
|
template<class M > |
int | Dune::countNonZeros (const M &matrix) |
| Get the number of nonzero fields in the matrix. More...
|
|
Matrix and Vector classes that support a block recursive structure capable of representing the natural structure from Finite Element discretisations.
The interface of our matrices is designed according to what they represent from a mathematical point of view. The vector classes are representations of vector spaces:
- FieldVector represents a vector space \(V=K^n\) where the field \(K\) is represented by a numeric type (e.g. double, float, complex). \(n\) is known at compile time.
- BlockVector represents a vector space \(V=W\times W \times W \times\cdots\times W\) where W is itself a vector space.
- VariableBlockVector represents a vector space having a two-level block structure of the form \(V=B^{n_1}\times B^{n_2}\times\ldots \times B^{n_m}\), i.e. it is constructed from \(m\) vector spaces, \(i=1,\ldots,m\).
The matrix classes represent linear maps \(A: V \mapsto W\) from vector space \(V\) to vector space \(W\) the recursive block structure of the matrix rows and columns immediately follows from the recursive block structure of the vectors representing the domain and range of the mapping, respectively:
- FieldMatrix represents a linear map \(M: V_1 \to V_2\) where \(V_1=K^n\) and \(V_2=K^m\) are vector spaces over the same field represented by a numerix type.
- BCRSMatrix represents a linear map \(M: V_1 \to V_2\) where \(V_1=W\times W \times W \times\cdots\times W\) and \(V_2=W\times W \times W \times\cdots\times W\) where W is itself a vector space.
- VariableBCRSMatrix is not yet implemented.
◆ countNonZeros()
template<class M >
int Dune::countNonZeros |
( |
const M & |
matrix | ) |
|
|
inline |
Get the number of nonzero fields in the matrix.
This is not the number of nonzero blocks, but the number of non zero scalar entries (on blocklevel 1) if the matrix is viewed as a flat matrix.
For FieldMatrix this is simply the number of columns times the number of rows, for a BCRSMatrix<FieldMatrix<K,n,m>> this is the number of nonzero blocks time n*m.
Referenced by Dune::Amg::MatrixHierarchy< M, PI, A >::build().
◆ matMultMat()
template<class T , class A , class A1 , class A2 , int n, int m, int k>
void Dune::matMultMat |
( |
BCRSMatrix< FieldMatrix< T, n, m >, A > & |
res, |
|
|
const BCRSMatrix< FieldMatrix< T, n, k >, A1 > & |
mat, |
|
|
const BCRSMatrix< FieldMatrix< T, k, m >, A2 > & |
matt, |
|
|
bool |
tryHard = false |
|
) |
| |
Calculate product of two sparse matrices ( \(C=A*B\)).
- Parameters
-
◆ matMultTransposeMat()
template<class T , class A , class A1 , class A2 , int n, int m, int k>
void Dune::matMultTransposeMat |
( |
BCRSMatrix< FieldMatrix< T, n, k >, A > & |
res, |
|
|
const BCRSMatrix< FieldMatrix< T, n, m >, A1 > & |
mat, |
|
|
const BCRSMatrix< FieldMatrix< T, k, m >, A2 > & |
matt, |
|
|
bool |
tryHard = false |
|
) |
| |
Calculate product of a sparse matrix with a transposed sparse matrices ( \(C=A*B^T\)).
- Parameters
-
res | Matrix for the result of the computation. |
mat | Matrix A. |
matt | Matrix B, which will be transposed before the multiplication. |
tryHard | ignored |
◆ transposeMatMultMat()
template<class T , class A , class A1 , class A2 , int n, int m, int k>
void Dune::transposeMatMultMat |
( |
BCRSMatrix< FieldMatrix< T, n, m >, A > & |
res, |
|
|
const BCRSMatrix< FieldMatrix< T, k, n >, A1 > & |
mat, |
|
|
const BCRSMatrix< FieldMatrix< T, k, m >, A2 > & |
matt, |
|
|
bool |
tryHard = false |
|
) |
| |
Calculate product of a transposed sparse matrix with another sparse matrices ( \(C=A^T*B\)).
- Parameters
-
res | Matrix for the result of the computation. |
mat | Matrix A, which will be transposed before the multiplication. |
matt | Matrix B. |
tryHard | ignored |