Dune Core Modules (2.9.0)
Matrix and Vector classes that support a block recursive structure capable of representing the natural structure from Finite Element discretisations. More...
Modules | |
Block Recursive Iterative Kernels | |
IO for matrices and vectors. | |
Provides methods for reading and writing matrices and vectors in various formats. | |
Dense Matrix and Vector Template Library | |
Type traits to retrieve the field and the real type of classes. | |
Files | |
file | matrixmatrix.hh |
provides functions for sparse matrix matrix multiplication. | |
file | matrixutils.hh |
Some handy generic functions for ISTL matrices. | |
Classes | |
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::MultiTypeBlockMatrix< FirstRow, Args > |
A Matrix class to support different block types. More... | |
class | Dune::MultiTypeBlockMatrix_Solver_Col< I, crow, ccol, remain_col > |
part of solvers for MultiTypeBlockVector & MultiTypeBlockMatrix types More... | |
class | Dune::MultiTypeBlockMatrix_Solver< I, crow, remain_row > |
solver for MultiTypeBlockVector & MultiTypeBlockMatrix types More... | |
struct | std::tuple_element< i, Dune::MultiTypeBlockMatrix< Args... > > |
Make std::tuple_element work for MultiTypeBlockMatrix. More... | |
class | Dune::MultiTypeBlockVector< Args > |
A Vector class to support different block types. More... | |
struct | std::tuple_element< i, Dune::MultiTypeBlockVector< Args... > > |
Make std::tuple_element work for MultiTypeBlockVector. More... | |
class | Dune::VariableBlockVector< B, A > |
A Vector of blocks with different blocksizes. More... | |
Typedefs | |
typedef MultiTypeBlockMatrix< FirstRow, Args... > | Dune::MultiTypeBlockMatrix< FirstRow, Args >::type |
using | Dune::MultiTypeBlockMatrix< FirstRow, Args >::size_type = std::size_t |
Type used for sizes. | |
using | Dune::MultiTypeBlockVector< Args >::size_type = std::size_t |
Type used for vector sizes. | |
typedef MultiTypeBlockVector< Args... > | Dune::MultiTypeBlockVector< Args >::type |
using | Dune::MultiTypeBlockVector< Args >::field_type = Std::detected_t< std::common_type_t, typename FieldTraits< std::decay_t< Args > >::field_type... > |
The type used for scalars. More... | |
Functions | |
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 > | |
auto | Dune::countNonZeros (const M &, typename std::enable_if_t< Dune::IsNumber< M >::value > *sfinae=nullptr) |
Get the number of nonzero fields in the matrix. More... | |
static constexpr size_type | Dune::MultiTypeBlockMatrix< FirstRow, Args >::N () |
Return the number of matrix rows. | |
static constexpr size_type | Dune::MultiTypeBlockMatrix< FirstRow, Args >::size () |
Return the number of matrix rows. More... | |
static constexpr size_type | Dune::MultiTypeBlockMatrix< FirstRow, Args >::M () |
Return the number of matrix columns. | |
template<size_type index> | |
auto | Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator[] (const std::integral_constant< size_type, index > indexVariable) -> decltype(std::get< index >(*this)) |
Random-access operator. More... | |
template<size_type index> | |
auto | Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator[] (const std::integral_constant< size_type, index > indexVariable) const -> decltype(std::get< index >(*this)) |
Const random-access operator. More... | |
template<typename T > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator= (const T &newval) |
MultiTypeBlockMatrix & | Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator*= (const field_type &k) |
vector space multiplication with scalar | |
MultiTypeBlockMatrix & | Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator/= (const field_type &k) |
vector space division by scalar | |
MultiTypeBlockMatrix & | Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator+= (const MultiTypeBlockMatrix &b) |
Add the entries of another matrix to this one. More... | |
MultiTypeBlockMatrix & | Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator-= (const MultiTypeBlockMatrix &b) |
Subtract the entries of another matrix from this one. More... | |
template<typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::mv (const X &x, Y &y) const |
y = A x | |
template<typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::umv (const X &x, Y &y) const |
y += A x | |
template<typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::mmv (const X &x, Y &y) const |
y -= A x | |
template<typename AlphaType , typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::usmv (const AlphaType &alpha, const X &x, Y &y) const |
y += alpha A x | |
template<typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::mtv (const X &x, Y &y) const |
y = A^T x | |
template<typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::umtv (const X &x, Y &y) const |
y += A^T x | |
template<typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::mmtv (const X &x, Y &y) const |
y -= A^T x | |
template<typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::usmtv (const field_type &alpha, const X &x, Y &y) const |
y += alpha A^T x | |
template<typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::umhv (const X &x, Y &y) const |
y += A^H x | |
template<typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::mmhv (const X &x, Y &y) const |
y -= A^H x | |
template<typename X , typename Y > | |
void | Dune::MultiTypeBlockMatrix< FirstRow, Args >::usmhv (const field_type &alpha, const X &x, Y &y) const |
y += alpha A^H x | |
auto | Dune::MultiTypeBlockMatrix< FirstRow, Args >::frobenius_norm2 () const |
square of frobenius norm, need for block recursion | |
FieldTraits< field_type >::real_type | Dune::MultiTypeBlockMatrix< FirstRow, Args >::frobenius_norm () const |
frobenius norm: sqrt(sum over squared values of entries) | |
auto | Dune::MultiTypeBlockMatrix< FirstRow, Args >::infinity_norm () const |
Bastardized version of the infinity-norm / row-sum norm. | |
auto | Dune::MultiTypeBlockMatrix< FirstRow, Args >::infinity_norm_real () const |
Bastardized version of the infinity-norm / row-sum norm. | |
template<typename T1 , typename... Args> | |
std::ostream & | Dune::operator<< (std::ostream &s, const MultiTypeBlockMatrix< T1, Args... > &m) |
<< operator for a MultiTypeBlockMatrix More... | |
template<typename Trhs , typename TVector , typename TMatrix , typename K > | |
static void | Dune::MultiTypeBlockMatrix_Solver_Col< I, crow, ccol, remain_col >::calc_rhs (const TMatrix &A, TVector &x, TVector &v, Trhs &b, const K &w) |
template<typename TVector , typename TMatrix , typename K > | |
static void | Dune::MultiTypeBlockMatrix_Solver< I, crow, remain_row >::dbgs (const TMatrix &A, TVector &x, const TVector &b, const K &w) |
template<typename TVector , typename TMatrix , typename K > | |
static void | Dune::MultiTypeBlockMatrix_Solver< I, crow, remain_row >::bsorf (const TMatrix &A, TVector &x, const TVector &b, const K &w) |
template<typename TVector , typename TMatrix , typename K > | |
static void | Dune::MultiTypeBlockMatrix_Solver< I, crow, remain_row >::bsorb (const TMatrix &A, TVector &x, const TVector &b, const K &w) |
template<typename TVector , typename TMatrix , typename K > | |
static void | Dune::MultiTypeBlockMatrix_Solver< I, crow, remain_row >::dbjac (const TMatrix &A, TVector &x, const TVector &b, const K &w) |
static constexpr size_type | Dune::MultiTypeBlockVector< Args >::size () |
Return the number of non-zero vector entries. More... | |
static constexpr size_type | Dune::MultiTypeBlockVector< Args >::N () |
Number of elements. | |
int | Dune::MultiTypeBlockVector< Args >::count () const |
size_type | Dune::MultiTypeBlockVector< Args >::dim () const |
Number of scalar elements. | |
template<size_type index> | |
std::tuple_element< index, TupleType >::type & | Dune::MultiTypeBlockVector< Args >::operator[] (const std::integral_constant< size_type, index > indexVariable) |
Random-access operator. More... | |
template<size_type index> | |
const std::tuple_element< index, TupleType >::type & | Dune::MultiTypeBlockVector< Args >::operator[] (const std::integral_constant< size_type, index > indexVariable) const |
Const random-access operator. More... | |
template<typename T > | |
void | Dune::MultiTypeBlockVector< Args >::operator= (const T &newval) |
Assignment operator. | |
void | Dune::MultiTypeBlockVector< Args >::operator+= (const type &newv) |
void | Dune::MultiTypeBlockVector< Args >::operator-= (const type &newv) |
template<class T , std::enable_if_t< IsNumber< T >::value, int > = 0> | |
void | Dune::MultiTypeBlockVector< Args >::operator*= (const T &w) |
Multiplication with a scalar. | |
template<class T , std::enable_if_t< IsNumber< T >::value, int > = 0> | |
void | Dune::MultiTypeBlockVector< Args >::operator/= (const T &w) |
Division by a scalar. | |
auto | Dune::MultiTypeBlockVector< Args >::one_norm () const |
Compute the 1-norm. | |
auto | Dune::MultiTypeBlockVector< Args >::one_norm_real () const |
Compute the simplified 1-norm (uses 1-norm also for complex values) | |
FieldTraits< field_type >::real_type | Dune::MultiTypeBlockVector< Args >::two_norm2 () const |
Compute the squared Euclidean norm. | |
FieldTraits< field_type >::real_type | Dune::MultiTypeBlockVector< Args >::two_norm () const |
Compute the Euclidean norm. | |
FieldTraits< field_type >::real_type | Dune::MultiTypeBlockVector< Args >::infinity_norm () const |
Compute the maximum norm. | |
FieldTraits< field_type >::real_type | Dune::MultiTypeBlockVector< Args >::infinity_norm_real () const |
Compute the simplified maximum norm (uses 1-norm for complex values) | |
template<typename Ta > | |
void | Dune::MultiTypeBlockVector< Args >::axpy (const Ta &a, const type &y) |
Axpy operation on this vector (*this += a * y) More... | |
template<typename... Args> | |
std::ostream & | Dune::operator<< (std::ostream &s, const MultiTypeBlockVector< Args... > &v) |
Send MultiTypeBlockVector to an outstream. | |
Detailed Description
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.
Typedef Documentation
◆ field_type
using Dune::MultiTypeBlockVector< Args >::field_type = Std::detected_t<std::common_type_t, typename FieldTraits< std::decay_t<Args> >::field_type...> |
The type used for scalars.
Use the std::common_type_t
of the Args
' field_type and use nonesuch
if no implementation of std::common_type
is provided for the given field_type
arguments.
◆ type [1/2]
typedef MultiTypeBlockMatrix<FirstRow, Args...> Dune::MultiTypeBlockMatrix< FirstRow, Args >::type |
own class' type
◆ type [2/2]
typedef MultiTypeBlockVector<Args...> Dune::MultiTypeBlockVector< Args >::type |
own class' type
Function Documentation
◆ axpy()
|
inline |
Axpy operation on this vector (*this += a * y)
- Template Parameters
-
Ta Type of the scalar 'a'
References Dune::Hybrid::forEach(), and Dune::Hybrid::integralRange().
◆ bsorb()
|
inlinestatic |
bsorb for MultiTypeBlockMatrix & MultiTypeBlockVector types
◆ bsorf()
|
inlinestatic |
bsorf for MultiTypeBlockMatrix & MultiTypeBlockVector types
◆ calc_rhs()
|
inlinestatic |
iterating over one row in MultiTypeBlockMatrix to calculate right side b-A[i][j]*x[j]
◆ count()
|
inline |
number of elements
- Deprecated:
- Use method
N
instead. This will be removed after Dune 2.8.
◆ countNonZeros()
|
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().
◆ dbgs()
|
inlinestatic |
Gauss-Seidel solver for MultiTypeBlockMatrix & MultiTypeBlockVector types
◆ dbjac()
|
inlinestatic |
Jacobi solver for MultiTypeBlockMatrix & MultiTypeBlockVector types
◆ matMultMat()
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 |
||
) |
◆ matMultTransposeMat()
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 |
||
) |
◆ operator+=() [1/2]
|
inline |
Add the entries of another matrix to this one.
- Parameters
-
b The matrix to add to this one. Its sparsity pattern has to be subset of the sparsity pattern of this matrix.
References Dune::Hybrid::forEach(), Dune::Hybrid::integralRange(), and Dune::MultiTypeBlockMatrix< FirstRow, Args >::size().
◆ operator+=() [2/2]
|
inline |
operator for MultiTypeBlockVector += MultiTypeBlockVector operations
References Dune::Hybrid::forEach(), and Dune::Hybrid::integralRange().
◆ operator-=() [1/2]
|
inline |
Subtract the entries of another matrix from this one.
- Parameters
-
b The matrix to subtract from this one. Its sparsity pattern has to be subset of the sparsity pattern of this matrix.
References Dune::Hybrid::forEach(), Dune::Hybrid::integralRange(), and Dune::MultiTypeBlockMatrix< FirstRow, Args >::size().
◆ operator-=() [2/2]
|
inline |
operator for MultiTypeBlockVector -= MultiTypeBlockVector operations
References Dune::Hybrid::forEach(), and Dune::Hybrid::integralRange().
◆ operator<<()
std::ostream & Dune::operator<< | ( | std::ostream & | s, |
const MultiTypeBlockMatrix< T1, Args... > & | m | ||
) |
<< operator for a MultiTypeBlockMatrix
operator<< for printing out a MultiTypeBlockMatrix
References Dune::Hybrid::forEach(), and Dune::Hybrid::integralRange().
◆ operator=()
|
inline |
assignment operator
References Dune::Hybrid::forEach(), Dune::Hybrid::integralRange(), and Dune::MultiTypeBlockMatrix< FirstRow, Args >::size().
◆ operator[]() [1/4]
|
inline |
Random-access operator.
This method mimicks the behavior of normal vector access with square brackets like, e.g., v[5] = 1. The problem is that the return type is different for each value of the argument in the brackets. Therefore we implement a trick using std::integral_constant. To access the first entry of a MultiTypeBlockVector named v write
The name '_0' used here as a static replacement of the integer number zero is arbitrary. Any other variable name can be used. If you don't like the separate variable, you can writee
◆ operator[]() [2/4]
|
inline |
Random-access operator.
This method mimicks the behavior of normal vector access with square brackets like, e.g., m[5] = .... The problem is that the return type is different for each value of the argument in the brackets. Therefore we implement a trick using std::integral_constant. To access the first row of a MultiTypeBlockMatrix named m write
The name '_0' used here as a static replacement of the integer number zero is arbitrary. Any other variable name can be used. If you don't like the separate variable, you can write
◆ operator[]() [3/4]
|
inline |
Const random-access operator.
This is the const version of the random-access operator. See the non-const version for a full explanation of how to use it.
◆ operator[]() [4/4]
|
inline |
Const random-access operator.
This is the const version of the random-access operator. See the non-const version for a full explanation of how to use it.
◆ size() [1/2]
|
inlinestaticconstexpr |
Return the number of matrix rows.
- Deprecated:
- Use method
N
instead. This will be removed after Dune 2.8.
Referenced by Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator*=(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator+=(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator-=(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator/=(), and Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator=().
◆ size() [2/2]
|
inlinestaticconstexpr |
Return the number of non-zero vector entries.
As this is a dense vector data structure, all entries are non-zero, and hence 'size' returns the same number as 'N'.
◆ transposeMatMultMat()
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
Referenced by MatrixInfo< BCRSMatrix >::computeNonSymCond2().