dune-istl 2.10
Loading...
Searching...
No Matches
Sparse Matrix and Vector classes

Matrix and Vector classes that support a block recursive structure capable of representing the natural structure from Finite Element discretisations. More...

Collaboration diagram for Sparse Matrix and Vector classes:

Topics

 Block Recursive Iterative Kernels
 IO for matrices and vectors.
 Provides methods for reading and writing matrices and vectors in various formats.
 DenseMatVec

Files

file  matrixmatrix.hh
 provides functions for sparse matrix matrix multiplication.
file  matrixutils.hh
 Some handy generic functions for ISTL matrices.

Classes

struct  Dune::MatrixDimension< M >
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::ImplicitMatrixBuilder< M_ >::row_object
 Proxy row object for entry access. More...
class  Dune::RealRowIterator< T >
class  Dune::BDMatrix< B, A >
 A block-diagonal matrix. More...
struct  Dune::FieldTraits< BDMatrix< B, A > >
class  Dune::BTDMatrix< B, A >
 A block-tridiagonal matrix. More...
struct  Dune::FieldTraits< BTDMatrix< B, A > >
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::FieldTraits< Matrix< T, A > >
struct  Dune::MatMultMatResult< M1, M2 >
 Helper TMP to get the result type of a sparse matrix matrix multiplication ( $C=A*B$). More...
struct  Dune::MatMultMatResult< FieldMatrix< T, n, k >, FieldMatrix< T, k, m > >
struct  Dune::MatMultMatResult< BCRSMatrix< FieldMatrix< T, n, k >, A >, BCRSMatrix< FieldMatrix< T, k, m >, A1 > >
struct  Dune::TransposedMatMultMatResult< M1, M2 >
 Helper TMP to get the result type of a sparse matrix matrix multiplication ( $C=A*B$). More...
struct  Dune::TransposedMatMultMatResult< FieldMatrix< T, k, n >, FieldMatrix< T, k, m > >
struct  Dune::TransposedMatMultMatResult< BCRSMatrix< FieldMatrix< T, k, n >, A >, BCRSMatrix< FieldMatrix< T, k, m >, A1 > >
struct  Dune::CheckIfDiagonalPresent< Matrix, blocklevel, l >
 Check whether the a matrix has diagonal values on blocklevel recursion levels. More...
struct  Dune::CheckIfDiagonalPresent< Matrix, 0, l >
class  Dune::MultiTypeBlockMatrix< FirstRow, Args >
 A Matrix class to support different block types. More...
struct  Dune::CheckIfDiagonalPresent< MultiTypeBlockMatrix< T1, Args... >, blocklevel, l >
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...
struct  std::tuple_size< Dune::MultiTypeBlockVector< Args... > >
 Make std::tuple_size work for MultiTypeBlockVector. More...
class  Dune::VariableBlockVector< B, A >
 A Vector of blocks with different blocksizes. More...

Typedefs

typedef M_ Dune::ImplicitMatrixBuilder< M_ >::Matrix
 The underlying matrix.
typedef Matrix::block_type Dune::ImplicitMatrixBuilder< M_ >::block_type
 The block_type of the underlying matrix.
typedef Matrix::size_type Dune::ImplicitMatrixBuilder< M_ >::size_type
 The size_type of the underlying matrix.
using Dune::field_type = typename Imp::BlockTraits<B>::field_type
typedef BCRSMatrix< FieldMatrix< T, n, m >, A >::CreateIterator Dune::SparsityPatternInitializer< T, A, n, m >::CreateIterator
typedef BCRSMatrix< FieldMatrix< T, n, m >, A >::size_type Dune::SparsityPatternInitializer< T, A, n, m >::size_type
typedef Dune::BCRSMatrix< FieldMatrix< T, n, m >, TA > Dune::MatrixInitializer< transpose, T, TA, n, m >::Matrix
typedef Matrix::CreateIterator Dune::MatrixInitializer< transpose, T, TA, n, m >::CreateIterator
typedef Matrix::size_type Dune::MatrixInitializer< transpose, T, TA, n, m >::size_type
typedef Dune::BCRSMatrix< Dune::FieldMatrix< T, n, m >, TA > Dune::MatrixInitializer< 1, T, TA, n, m >::Matrix
typedef Matrix::CreateIterator Dune::MatrixInitializer< 1, T, TA, n, m >::CreateIterator
typedef Matrix::size_type Dune::MatrixInitializer< 1, T, TA, n, m >::size_type
typedef BCRSMatrix< FieldMatrix< T, n, m >, A > Dune::EntryAccumulatorFather< T, A, n, m >::Matrix
typedef Matrix::RowIterator Dune::EntryAccumulatorFather< T, A, n, m >::Row
typedef Matrix::ColIterator Dune::EntryAccumulatorFather< T, A, n, m >::Col
typedef BCRSMatrix< FieldMatrix< T, n, m >, A > Dune::EntryAccumulator< T, A, n, m, transpose >::Matrix
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, transpose >::size_type
typedef BCRSMatrix< FieldMatrix< T, n, m >, A > Dune::EntryAccumulator< T, A, n, m, 0 >::Matrix
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 0 >::size_type
typedef BCRSMatrix< FieldMatrix< T, n, m >, A > Dune::EntryAccumulator< T, A, n, m, 1 >::Matrix
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 1 >::size_type
typedef BCRSMatrix< FieldMatrix< T, n, m >, A > Dune::EntryAccumulator< T, A, n, m, 2 >::Matrix
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 2 >::size_type
typedef FieldMatrix< T, n, m > Dune::MatMultMatResult< FieldMatrix< T, n, k >, FieldMatrix< T, k, m > >::type
typedef BCRSMatrix< typename MatMultMatResult< FieldMatrix< T, n, k >, FieldMatrix< T, k, m > >::type, std::allocator< typename MatMultMatResult< FieldMatrix< T, n, k >, FieldMatrix< T, k, m > >::type > > Dune::MatMultMatResult< BCRSMatrix< FieldMatrix< T, n, k >, A >, BCRSMatrix< FieldMatrix< T, k, m >, A1 > >::type
typedef FieldMatrix< T, n, m > Dune::TransposedMatMultMatResult< FieldMatrix< T, k, n >, FieldMatrix< T, k, m > >::type
typedef BCRSMatrix< typename MatMultMatResult< FieldMatrix< T, n, k >, FieldMatrix< T, k, m > >::type, std::allocator< typename MatMultMatResult< FieldMatrix< T, n, k >, FieldMatrix< T, k, m > >::type > > Dune::TransposedMatMultMatResult< BCRSMatrix< FieldMatrix< T, k, n >, A >, BCRSMatrix< FieldMatrix< T, k, m >, A1 > >::type
using Dune::FieldTraits< MultiTypeBlockVector< Args... > >::field_type = typename MultiTypeBlockVector<Args...>::field_type
using Dune::FieldTraits< MultiTypeBlockVector< Args... > >::real_type = typename MultiTypeBlockVector<Args...>::real_type
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.
using Dune::MultiTypeBlockVector< Args >::real_type = Std::detected_t<std::common_type_t, typename FieldTraits< std::decay_t<Args> >::real_type...>
 The type used for real values.
using std::tuple_element< i, Dune::MultiTypeBlockVector< Args... > >::type = typename std::tuple_element<i, std::tuple<Args...> >::type

Enumerations

enum  Dune::BuildMode { Dune::row_wise , Dune::random , Dune::implicit , Dune::unknown }
enum  { Dune::SparsityPatternInitializer< T, A, n, m >::do_break =true }
enum  { Dune::MatrixInitializer< transpose, T, TA, n, m >::do_break =true }
enum  { Dune::MatrixInitializer< 1, T, TA, n, m >::do_break =false }
enum  { Dune::EntryAccumulatorFather< T, A, n, m >::do_break =false }

Functions

block_typeDune::ImplicitMatrixBuilder< M_ >::row_object::operator[] (size_type j) const
 Returns entry in column j.
 Dune::ImplicitMatrixBuilder< M_ >::ImplicitMatrixBuilder (Matrix &m)
 Creates an ImplicitMatrixBuilder for matrix m.
 Dune::ImplicitMatrixBuilder< M_ >::ImplicitMatrixBuilder (Matrix &m, size_type rows, size_type cols, size_type avg_cols_per_row, double overflow_fraction)
 Sets up matrix m for implicit construction using the given parameters and creates an ImplicitBmatrixuilder for it.
row_object Dune::ImplicitMatrixBuilder< M_ >::operator[] (size_type i) const
 Returns a proxy for entries in row i.
size_type Dune::ImplicitMatrixBuilder< M_ >::N () const
 The number of rows in the matrix.
size_type Dune::ImplicitMatrixBuilder< M_ >::M () const
 The number of columns in the matrix.
*random access to the rows row_typeDune::operator[] (size_type i)
*constructor Dune::RealRowIterator< T >::RealRowIterator (row_type *_p, size_type _i)
*empty use with care Dune::RealRowIterator< T >::!RealRowIterator ()
 Dune::RealRowIterator< T >::RealRowIterator (const RealRowIterator< ValueType > &it)
*return index size_type Dune::RealRowIterator< T >::index () const
std::ptrdiff_t Dune::RealRowIterator< T >::distanceTo (const RealRowIterator< ValueType > &other) const
std::ptrdiff_t Dune::RealRowIterator< T >::distanceTo (const RealRowIterator< const ValueType > &other) const
*equality bool Dune::RealRowIterator< T >::equals (const RealRowIterator< ValueType > &other) const
*equality bool Dune::RealRowIterator< T >::equals (const RealRowIterator< const ValueType > &other) const
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$).
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$).
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$).
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.
template<class M>
auto Dune::countNonZeros (const M &matrix, typename std::enable_if_t<!Dune::IsNumber< M >::value > *sfinae=nullptr)
static constexpr size_type Dune::MultiTypeBlockVector< Args >::size ()
 Return the number of non-zero vector entries.
static constexpr size_type Dune::MultiTypeBlockVector< Args >::N ()
 Number of elements.
size_type Dune::MultiTypeBlockVector< Args >::dim () const
 Number of scalar elements.
template<size_type index>
std::tuple_element< index, TupleType >::typeDune::MultiTypeBlockVector< Args >::operator[] (const std::integral_constant< size_type, index > indexVariable)
 Random-access operator.
template<size_type index>
const std::tuple_element< index, TupleType >::typeDune::MultiTypeBlockVector< Args >::operator[] (const std::integral_constant< size_type, index > indexVariable) const
 Const random-access operator.
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.
field_type Dune::MultiTypeBlockVector< Args >::operator* (const type &newv) const
field_type Dune::MultiTypeBlockVector< Args >::dot (const type &newv) const
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).
real_type Dune::MultiTypeBlockVector< Args >::two_norm2 () const
 Compute the squared Euclidean norm.
real_type Dune::MultiTypeBlockVector< Args >::two_norm () const
 Compute the Euclidean norm.
real_type Dune::MultiTypeBlockVector< Args >::infinity_norm () const
 Compute the maximum norm.
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).
template<typename... Args>
std::ostream & Dune::operator<< (std::ostream &s, const MultiTypeBlockVector< Args... > &v)
 Send MultiTypeBlockVector to an outstream.

Variables

*the type representing the components typedef B Dune::block_type
*the allocator type typedef A Dune::allocator_type
*The type for the index access and the size typedef A::size_type Dune::size_type
*implement row_type with compressed vector typedef Imp::CompressedBlockVectorWindow< B, size_typeDune::row_type
*brief The unqualified value type typedef std::remove_const< T >::type Dune::RealRowIterator< T >::ValueType
double Dune::CompressionStatistics< size_type >::avg
 average number of non-zeroes per row.
size_type Dune::CompressionStatistics< size_type >::maximum
 maximum number of non-zeroes per row.
size_type Dune::CompressionStatistics< size_type >::overflow_total
 total number of elements written to the overflow area during construction.
double Dune::CompressionStatistics< size_type >::mem_ratio
 fraction of wasted memory resulting from non-used overflow area.
 Dune::notbuilt =0
 A sparse block matrix with compressed row storage.
 Dune::notAllocated =0
 Matrix is not built at all, no memory has been allocated, build mode and size can still be set.
 Dune::building =1
 Matrix is currently being built, some memory has been allocated, build mode and size are fixed.
 Dune::rowSizesBuilt =2
 The row sizes of the matrix are known.
 Dune::built
 The matrix structure is fully built.
*empty Dune::RealRowIterator< T >::constructor
MatrixDune::EntryAccumulatorFather< T, A, n, m >::mat
Col Dune::EntryAccumulatorFather< T, A, n, m >::col

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:

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:

Typedef Documentation

◆ block_type

template<class M_>
typedef Matrix::block_type Dune::ImplicitMatrixBuilder< M_ >::block_type

The block_type of the underlying matrix.

◆ Col

template<class T, class A, int n, int m>
typedef Matrix::ColIterator Dune::EntryAccumulatorFather< T, A, n, m >::Col

◆ CreateIterator [1/3]

template<int transpose, class T, class TA, int n, int m>
typedef Matrix::CreateIterator Dune::MatrixInitializer< transpose, T, TA, n, m >::CreateIterator

◆ CreateIterator [2/3]

template<class T, class TA, int n, int m>
typedef Matrix::CreateIterator Dune::MatrixInitializer< 1, T, TA, n, m >::CreateIterator

◆ CreateIterator [3/3]

template<class T, class A, int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A>::CreateIterator Dune::SparsityPatternInitializer< T, A, n, m >::CreateIterator

◆ field_type [1/3]

using Dune::field_type = typename Imp::BlockTraits<B>::field_type

◆ field_type [2/3]

template<typename... Args>
using Dune::FieldTraits< MultiTypeBlockVector< Args... > >::field_type = typename MultiTypeBlockVector<Args...>::field_type

◆ field_type [3/3]

template<typename... Args>
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.

◆ Matrix [1/8]

template<class T, class A, int n, int m, int transpose>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Dune::EntryAccumulator< T, A, n, m, transpose >::Matrix

◆ Matrix [2/8]

template<class T, class A, int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Dune::EntryAccumulator< T, A, n, m, 0 >::Matrix

◆ Matrix [3/8]

template<class T, class A, int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Dune::EntryAccumulator< T, A, n, m, 1 >::Matrix

◆ Matrix [4/8]

template<class T, class A, int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Dune::EntryAccumulator< T, A, n, m, 2 >::Matrix

◆ Matrix [5/8]

template<class T, class A, int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Dune::EntryAccumulatorFather< T, A, n, m >::Matrix

◆ Matrix [6/8]

template<class M_>
typedef M_ Dune::ImplicitMatrixBuilder< M_ >::Matrix

The underlying matrix.

◆ Matrix [7/8]

template<int transpose, class T, class TA, int n, int m>
typedef Dune::BCRSMatrix<FieldMatrix<T,n,m>,TA> Dune::MatrixInitializer< transpose, T, TA, n, m >::Matrix

◆ Matrix [8/8]

template<class T, class TA, int n, int m>
typedef Dune::BCRSMatrix<Dune::FieldMatrix<T,n,m>,TA> Dune::MatrixInitializer< 1, T, TA, n, m >::Matrix

◆ real_type [1/2]

template<typename... Args>
using Dune::FieldTraits< MultiTypeBlockVector< Args... > >::real_type = typename MultiTypeBlockVector<Args...>::real_type

◆ real_type [2/2]

template<typename... Args>
using Dune::MultiTypeBlockVector< Args >::real_type = Std::detected_t<std::common_type_t, typename FieldTraits< std::decay_t<Args> >::real_type...>

The type used for real values.

Use the std::common_type_t of the Args' real_type and use nonesuch if no implementation of std::common_type is provided for the given real_type arguments.

◆ Row

template<class T, class A, int n, int m>
typedef Matrix::RowIterator Dune::EntryAccumulatorFather< T, A, n, m >::Row

◆ size_type [1/9]

template<class T, class A, int n, int m, int transpose>
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, transpose >::size_type

◆ size_type [2/9]

template<class T, class A, int n, int m>
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 0 >::size_type

◆ size_type [3/9]

template<class T, class A, int n, int m>
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 1 >::size_type

◆ size_type [4/9]

template<class T, class A, int n, int m>
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 2 >::size_type

◆ size_type [5/9]

template<class M_>
typedef Matrix::size_type Dune::ImplicitMatrixBuilder< M_ >::size_type

The size_type of the underlying matrix.

◆ size_type [6/9]

template<int transpose, class T, class TA, int n, int m>
typedef Matrix::size_type Dune::MatrixInitializer< transpose, T, TA, n, m >::size_type

◆ size_type [7/9]

template<class T, class TA, int n, int m>
typedef Matrix::size_type Dune::MatrixInitializer< 1, T, TA, n, m >::size_type

◆ size_type [8/9]

template<typename... Args>
using Dune::MultiTypeBlockVector< Args >::size_type = std::size_t

Type used for vector sizes.

◆ size_type [9/9]

template<class T, class A, int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A>::size_type Dune::SparsityPatternInitializer< T, A, n, m >::size_type

◆ type [1/6]

template<typename T, typename A, typename A1, int n, int k, int m>
typedef BCRSMatrix<typename MatMultMatResult<FieldMatrix<T,n,k>,FieldMatrix<T,k,m> >::type, std::allocator<typename MatMultMatResult<FieldMatrix<T,n,k>,FieldMatrix<T,k,m> >::type> > Dune::MatMultMatResult< BCRSMatrix< FieldMatrix< T, n, k >, A >, BCRSMatrix< FieldMatrix< T, k, m >, A1 > >::type

◆ type [2/6]

template<typename T, int n, int k, int m>
typedef FieldMatrix<T,n,m> Dune::MatMultMatResult< FieldMatrix< T, n, k >, FieldMatrix< T, k, m > >::type

◆ type [3/6]

template<typename... Args>
typedef MultiTypeBlockVector<Args...> Dune::MultiTypeBlockVector< Args >::type

own class' type

◆ type [4/6]

template<typename T, typename A, typename A1, int n, int k, int m>
typedef BCRSMatrix<typename MatMultMatResult<FieldMatrix<T,n,k>,FieldMatrix<T,k,m> >::type, std::allocator<typename MatMultMatResult<FieldMatrix<T,n,k>,FieldMatrix<T,k,m> >::type> > Dune::TransposedMatMultMatResult< BCRSMatrix< FieldMatrix< T, k, n >, A >, BCRSMatrix< FieldMatrix< T, k, m >, A1 > >::type

◆ type [5/6]

template<typename T, int n, int k, int m>
typedef FieldMatrix<T,n,m> Dune::TransposedMatMultMatResult< FieldMatrix< T, k, n >, FieldMatrix< T, k, m > >::type

◆ type [6/6]

template<size_t i, typename... Args>
using std::tuple_element< i, Dune::MultiTypeBlockVector< Args... > >::type = typename std::tuple_element<i, std::tuple<Args...> >::type

Enumeration Type Documentation

◆ anonymous enum

template<int transpose, class T, class TA, int n, int m>
anonymous enum
Enumerator
do_break 

◆ anonymous enum

template<class T, class TA, int n, int m>
anonymous enum
Enumerator
do_break 

◆ anonymous enum

template<class T, class A, int n, int m>
anonymous enum
Enumerator
do_break 

◆ anonymous enum

template<class T, class A, int n, int m>
anonymous enum
Enumerator
do_break 

◆ BuildMode

Enumerator
row_wise 

Build in a row-wise manner.

Rows are built up in sequential order. Size of the row and the column indices are defined. A row can be used as soon as it is initialized. With respect to memory there are two variants of this scheme: (a) number of non-zeroes known in advance (application finite difference schemes), (b) number of non-zeroes not known in advance (application: Sparse LU, ILU(n)).

random 

Build entries randomly.

For general finite element implementations the number of rows n is known, the number of non-zeroes might also be known (e.g. #edges + #nodes for P2) but the size of a row and the indices of a row cannot be defined in sequential order.

implicit 

Build entries randomly with an educated guess for the number of entries per row.

Allows random order generation as in random mode, but row sizes do not need to be given first. Instead an average number of non-zeroes per row is passed to the constructor. Matrix setup is finished with compress(), full data access during build stage is possible.

unknown 

Build mode not set!

Function Documentation

◆ !RealRowIterator()

template<class T>
*empty use with care Dune::RealRowIterator< T >::!RealRowIterator ( )
inline

◆ axpy()

template<typename... Args>
template<typename Ta>
void Dune::MultiTypeBlockVector< Args >::axpy ( const Ta & a,
const type & y )
inline

Axpy operation on this vector (*this += a * y).

Template Parameters
TaType of the scalar 'a'

◆ countNonZeros() [1/2]

template<class M>
auto Dune::countNonZeros ( const M & ,
typename std::enable_if_t< Dune::IsNumber< M >::value > * sfinae = nullptr )
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.

◆ countNonZeros() [2/2]

template<class M>
auto Dune::countNonZeros ( const M & matrix,
typename std::enable_if_t<!Dune::IsNumber< M >::value > * sfinae = nullptr )
inline

◆ dim()

template<typename... Args>
size_type Dune::MultiTypeBlockVector< Args >::dim ( ) const
inline

Number of scalar elements.

◆ distanceTo() [1/2]

template<class T>
std::ptrdiff_t Dune::RealRowIterator< T >::distanceTo ( const RealRowIterator< const ValueType > & other) const
inline

◆ distanceTo() [2/2]

template<class T>
std::ptrdiff_t Dune::RealRowIterator< T >::distanceTo ( const RealRowIterator< ValueType > & other) const
inline

◆ dot()

template<typename... Args>
field_type Dune::MultiTypeBlockVector< Args >::dot ( const type & newv) const
inline

◆ equals() [1/2]

template<class T>
*equality bool Dune::RealRowIterator< T >::equals ( const RealRowIterator< const ValueType > & other) const
inline

◆ equals() [2/2]

template<class T>
*equality bool Dune::RealRowIterator< T >::equals ( const RealRowIterator< ValueType > & other) const
inline

◆ ImplicitMatrixBuilder() [1/2]

template<class M_>
Dune::ImplicitMatrixBuilder< M_ >::ImplicitMatrixBuilder ( Matrix & m)
inline

Creates an ImplicitMatrixBuilder for matrix m.

Note
You can only pass a completely set up matrix to this constructor: All of setBuildMode(), setImplicitBuildModeParameters() and setSize() must have been called with the correct values.

◆ ImplicitMatrixBuilder() [2/2]

template<class M_>
Dune::ImplicitMatrixBuilder< M_ >::ImplicitMatrixBuilder ( Matrix & m,
size_type rows,
size_type cols,
size_type avg_cols_per_row,
double overflow_fraction )
inline

Sets up matrix m for implicit construction using the given parameters and creates an ImplicitBmatrixuilder for it.

Using this constructor, you can perform the necessary matrix setup and the creation of the ImplicitMatrixBuilder in a single step. The matrix must still be in the build stage notAllocated, otherwise a BCRSMatrixError will be thrown. For a detailed explanation of the matrix parameters, see BCRSMatrix.

Parameters
mthe matrix to be built
rowsthe number of matrix rows
colsthe number of matrix columns
avg_cols_per_rowthe average number of non-zero columns per row
overflow_fractionthe amount of overflow to reserve in the matrix
See also
BCRSMatrix

◆ index()

template<class T>
*return index size_type Dune::RealRowIterator< T >::index ( ) const
inline

◆ infinity_norm()

template<typename... Args>
real_type Dune::MultiTypeBlockVector< Args >::infinity_norm ( ) const
inline

Compute the maximum norm.

◆ infinity_norm_real()

template<typename... Args>
real_type Dune::MultiTypeBlockVector< Args >::infinity_norm_real ( ) const
inline

Compute the simplified maximum norm (uses 1-norm for complex values).

◆ M()

template<class M_>
size_type Dune::ImplicitMatrixBuilder< M_ >::M ( ) const
inline

The number of columns in the matrix.

◆ 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
resMatrix for the result of the computation.
matMatrix A.
mattMatrix B.
tryHardignored

◆ 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
resMatrix for the result of the computation.
matMatrix A.
mattMatrix B, which will be transposed before the multiplication.
tryHardignored

◆ N() [1/2]

template<class M_>
size_type Dune::ImplicitMatrixBuilder< M_ >::N ( ) const
inline

The number of rows in the matrix.

◆ N() [2/2]

template<typename... Args>
constexpr size_type Dune::MultiTypeBlockVector< Args >::N ( )
inlinestaticconstexpr

Number of elements.

◆ one_norm()

template<typename... Args>
auto Dune::MultiTypeBlockVector< Args >::one_norm ( ) const
inline

Compute the 1-norm.

◆ one_norm_real()

template<typename... Args>
auto Dune::MultiTypeBlockVector< Args >::one_norm_real ( ) const
inline

Compute the simplified 1-norm (uses 1-norm also for complex values).

◆ operator*()

template<typename... Args>
field_type Dune::MultiTypeBlockVector< Args >::operator* ( const type & newv) const
inline

◆ operator*=()

template<typename... Args>
template<class T, std::enable_if_t< IsNumber< T >::value, int > = 0>
void Dune::MultiTypeBlockVector< Args >::operator*= ( const T & w)
inline

Multiplication with a scalar.

◆ operator+=()

template<typename... Args>
void Dune::MultiTypeBlockVector< Args >::operator+= ( const type & newv)
inline

operator for MultiTypeBlockVector += MultiTypeBlockVector operations

◆ operator-=()

template<typename... Args>
void Dune::MultiTypeBlockVector< Args >::operator-= ( const type & newv)
inline

operator for MultiTypeBlockVector -= MultiTypeBlockVector operations

◆ operator/=()

template<typename... Args>
template<class T, std::enable_if_t< IsNumber< T >::value, int > = 0>
void Dune::MultiTypeBlockVector< Args >::operator/= ( const T & w)
inline

Division by a scalar.

◆ operator<<()

template<typename... Args>
std::ostream & Dune::operator<< ( std::ostream & s,
const MultiTypeBlockVector< Args... > & v )

Send MultiTypeBlockVector to an outstream.

◆ operator=()

template<typename... Args>
template<typename T>
void Dune::MultiTypeBlockVector< Args >::operator= ( const T & newval)
inline

Assignment operator.

◆ operator[]() [1/5]

template<class M_>
row_object Dune::ImplicitMatrixBuilder< M_ >::operator[] ( size_type i) const
inline

Returns a proxy for entries in row i.

◆ operator[]() [2/5]

template<class M_>
block_type & Dune::ImplicitMatrixBuilder< M_ >::row_object::operator[] ( size_type j) const
inline

Returns entry in column j.

◆ operator[]() [3/5]

template<typename... Args>
template<size_type index>
std::tuple_element< index, TupleType >::type & Dune::MultiTypeBlockVector< Args >::operator[] ( const std::integral_constant< size_type, index > indexVariable)
inline

Random-access operator.

This method mimics 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

std::integral_constant<std::size_t,0> _0;
v[_0] = ...
A Vector class to support different block types.
Definition multitypeblockvector.hh:59

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

v[std::integral_constant<std::size_t,0>()] = ...

◆ operator[]() [4/5]

template<typename... Args>
template<size_type index>
const std::tuple_element< index, TupleType >::type & Dune::MultiTypeBlockVector< Args >::operator[] ( const std::integral_constant< size_type, index > indexVariable) const
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[]() [5/5]

*same for read only access const row_type & Dune::operator[] ( size_type i)

◆ RealRowIterator() [1/2]

template<class T>
Dune::RealRowIterator< T >::RealRowIterator ( const RealRowIterator< ValueType > & it)
inline

◆ RealRowIterator() [2/2]

template<class T>
*constructor Dune::RealRowIterator< T >::RealRowIterator ( row_type * _p,
size_type _i )
inline

◆ size()

template<typename... Args>
constexpr size_type Dune::MultiTypeBlockVector< Args >::size ( )
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()

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
resMatrix for the result of the computation.
matMatrix A, which will be transposed before the multiplication.
mattMatrix B.
tryHardignored

◆ two_norm()

template<typename... Args>
real_type Dune::MultiTypeBlockVector< Args >::two_norm ( ) const
inline

Compute the Euclidean norm.

◆ two_norm2()

template<typename... Args>
real_type Dune::MultiTypeBlockVector< Args >::two_norm2 ( ) const
inline

Compute the squared Euclidean norm.

Variable Documentation

◆ allocator_type

* the allocator type typedef A Dune::allocator_type

◆ avg

template<typename size_type>
double Dune::CompressionStatistics< size_type >::avg

average number of non-zeroes per row.

◆ block_type

* the type representing the components typedef B Dune::block_type

◆ building

Dune::building =1

Matrix is currently being built, some memory has been allocated, build mode and size are fixed.

◆ built

Dune::built
Initial value:
=3
}

The matrix structure is fully built.

◆ col

template<class T, class A, int n, int m>
Col Dune::EntryAccumulatorFather< T, A, n, m >::col
protected

◆ constructor

template<class T>
* empty Dune::RealRowIterator< T >::constructor

◆ mat

template<class T, class A, int n, int m>
Matrix& Dune::EntryAccumulatorFather< T, A, n, m >::mat
protected

◆ maximum

template<typename size_type>
size_type Dune::CompressionStatistics< size_type >::maximum

maximum number of non-zeroes per row.

◆ mem_ratio

template<typename size_type>
double Dune::CompressionStatistics< size_type >::mem_ratio

fraction of wasted memory resulting from non-used overflow area.

mem_ratio is equal to nonzeros()/(# allocated matrix entries).

◆ notAllocated

Dune::notAllocated =0

Matrix is not built at all, no memory has been allocated, build mode and size can still be set.

◆ notbuilt

Dune::notbuilt =0

A sparse block matrix with compressed row storage.

Implements a block compressed row storage scheme. The block type B can be any type implementing the matrix interface.

Different ways to build up a compressed row storage matrix are supported:

  1. Row-wise scheme
  2. Random scheme
  3. implicit scheme

Error checking: no error checking is provided normally. Setting the compile time switch DUNE_ISTL_WITH_CHECKING enables error checking.

Details:

  1. Row-wise scheme

Rows are built up in sequential order. Size of the row and the column indices are defined. A row can be used as soon as it is initialized. With respect to memory there are two variants of this scheme: (a) number of non-zeroes known in advance (application finite difference schemes), (b) number of non-zeroes not known in advance (application: Sparse LU, ILU(n)).

#include<dune/common/fmatrix.hh>
...
// third parameter is an optional upper bound for the number
// of nonzeros. If given the matrix will use one array for all values
// as opposed to one for each row.
for(Iter row=B.createbegin(); row!=B.createend(); ++row){
// Add nonzeros for left neighbour, diagonal and right neighbour
if(row.index()>0)
row.insert(row.index()-1);
row.insert(row.index());
if(row.index()<B.N()-1)
row.insert(row.index()+1);
}
// Now the sparsity pattern is fully set up and we can add values
B[0][0]=2;
...
Implementation of the BCRSMatrix class.
Definition matrixutils.hh:24
Definition matrixutils.hh:27
  1. Random scheme

For general finite element implementations the number of rows n is known, the number of non-zeroes might also be known (e.g. #edges + #nodes for P2) but the size of a row and the indices of a row can not be defined in sequential order.

#include<dune/common/fmatrix.hh>
...
// initially set row size for each row
B.setrowsize(0,1);
B.setrowsize(3,4);
B.setrowsize(2,1);
B.setrowsize(1,1);
// increase row size for row 2
B.incrementrowsize(2);
// finalize row setup phase
B.endrowsizes();
// add column entries to rows
B.addindex(0,0);
B.addindex(3,1);
B.addindex(2,2);
B.addindex(1,1);
B.addindex(2,0);
B.addindex(3,2);
B.addindex(3,0);
B.addindex(3,3);
// finalize column setup phase
B.endindices();
// set entries using the random access operator
B[0][0] = 1;
B[1][1] = 2;
B[2][0] = 3;
B[2][2] = 4;
B[3][1] = 5;
B[3][2] = 6;
B[3][0] = 7;
B[3][3] = 8;
  1. implicit scheme

With the 'random scheme` described above, the sparsity pattern has to be determined and stored before the matrix is assembled. This requires a dedicated iteration over the grid elements, which can be costly in terms of time. Also, additional memory is needed to store the pattern before it can be given to the 'random' build mode.

On the other hand, often one has good a priori knowledge about the number of entries a row contains on average. The implicit mode tries to make use of that knowledge, and allows the setup of matrix pattern and numerical values together.

Constructing and filling a BCRSMatrix with the 'implicit' mode is performed in two steps: In a setup phase, matrix entries with numerical values can be inserted into the matrix. Then, a compression algorithm is called which defragments and optimizes the memory layout. After this compression step, the matrix is ready to be used, and no further nonzero entries can be added.

To use this mode, either construct a matrix object via

  • BCRSMatrix(size_type _n, size_type _m, size_type _avg, double compressionBufferSize, BuildMode bm)

or default-construct the matrix and then call

  • setImplicitBuildModeParameters(size_type _avg, double compressionBufferSize)
  • setSize(size_type rows, size_type columns, size_type nnz=0)

The parameter _avg specifies the expected number of (block) entries per matrix row.

When the BCRSMatrix object is first constructed with the 'implicit' build mode, two areas for matrix entry storage are allocated:

1) A large continuous chunk of memory that can hold the expected number of entries. In addition, this chunk contains an extra part of memory called the 'compression buffer', located before the memory for the matrix itself. The size of this buffer will be _avg * _n * compressionBufferSize.

2) An associative container indexed by $i,j$-pairs, which will hold surplus matrix entries during the setup phase (the 'overflow area'). Its content is merged into the main memory during compression.

You can then start filling your matrix by calling entry(size_type row, size_type col), which returns the corresponding matrix entry, creating it on the fly if it does not exist yet. The matrix pattern is hence created implicitly by simply accessing nonzero entries during the initial matrix assembly. Note that new entries are not zero-initialized, though, and hence the first operation on each entry has to be an assignment.

If a row contains more non-zero entries than what was specified in the _avg parameter, the surplus entries are stored in the 'overflow area' during the initial setup phase. After all indices are added, call compress() to trigger the compression step that optimizes the matrix and integrates any entries from the overflow area into the standard BCRS storage. This compression step builds up the final memory layout row by row. It will fail with an exception if the compression buffer is not large enough, which would lead to compressed rows overwriting uncompressed ones. More precisely, if $\textrm{nnz}_j$ denotes the number of non-zeros in the $j$-th row, then the compression algorithm will succeed if the maximal number of non-zeros in the $i$-th row is

\‍[      M_i = \textrm{avg} + A + \sum_{j<i} (\textrm{avg} - \textrm{nnz}_j)
\‍]

for all $i$, where $ A = \textrm{avg}(n \cdot \textrm{compressionBufferSize}) $ is the total size of the compression buffer determined by the parameters explained above.

The data of the matrix is now located at the beginning of the allocated area, and covers what used to be the compression buffer. In exchange, there is now unused space at the end of the large allocated piece of memory. This will go unused and cannot be freed during the lifetime of the matrix, but it has no negative impact on run-time performance. No matrix entries may be added after the compression step.

The compress() method returns a value of type Dune::CompressionStatistics, which you can inspect to tune the construction parameters _avg and compressionBufferSize.

Use of copy constructor, assignment operator and matrix vector arithmetic is not supported until the matrix is fully built.

The following sample code constructs a $ 10 \times 10$ matrix, with an expected number of two entries per matrix row. The compression buffer size is set to 0.4. Hence the main chunk of allocated memory will be able to hold 10 * 2 entries in the matrix rows, and 10 * 2 * 0.4 entries in the compression buffer. In total that's 28 entries.

M m(10, 10, 2, 0.4, M::implicit);
// Fill in some arbitrary entries; the order is irrelevant.
// Even operations on these would be possible, you get a reference to the entry!
m.entry(0,0) = 0.;
m.entry(8,0) = 0.;
m.entry(1,8) = 0.; m.entry(1,0) = 0.; m.entry(1,5) = 0.;
m.entry(2,0) = 0.;
m.entry(3,5) = 0.; m.entry(3,0) = 0.; m.entry(3,8) = 0.;
m.entry(4,0) = 0.;
m.entry(9,0) = 0.; m.entry(9,5) = 0.; m.entry(9,8) = 0.;
m.entry(5,0) = 0.; m.entry(5,5) = 0.; m.entry(5,8) = 0.;
m.entry(6,0) = 0.;
m.entry(7,0) = 0.; m.entry(7,5) = 0.; m.entry(7,8) = 0.;

Internally the index array now looks like this:

// xxxxxxxx0x800x500x050x050x05
// ........|.|.|.|.|.|.|.|.|.|.

The second row denotes the beginnings of the matrix rows. The eight 'x' on the left are the compression buffer. The overflow area contains the entries (1,5,0.0), (3,8,0.0), (5,8,0.0), (7,8,0.0), and (9,8,0.0). These are entries of rows 1, 3, 5, 7, and 9, which have three entries each, even though only two were anticipated.

//finish building by compressing the array
Statistics about compression achieved in implicit mode.
Definition bcrsmatrix.hh:88

Internally the index array now looks like this:

// 00580058005800580058xxxxxxxx
// ||..||..||..||..||..........

The compression buffer on the left is gone now, and the matrix has a real CRS layout. The 'x' on the right will be unused for the rest of the matrix' lifetime. */ template<class B, class A=std::allocator > class BCRSMatrix { friend struct MatrixDimension<BCRSMatrix>; public: enum BuildStage { /** Matrix is not built at all, no memory has been allocated, build mode and size can still be set.

◆ overflow_total

template<typename size_type>
size_type Dune::CompressionStatistics< size_type >::overflow_total

total number of elements written to the overflow area during construction.

◆ row_type

* implement row_type with compressed vector typedef Imp::CompressedBlockVectorWindow<B,size_type> Dune::row_type

◆ rowSizesBuilt

Dune::rowSizesBuilt =2

The row sizes of the matrix are known.

Only used in random mode.

◆ size_type

* The type for the index access and the size typedef A::size_type Dune::size_type

◆ ValueType

template<class T>
* brief The unqualified value type typedef std::remove_const<T>::type Dune::RealRowIterator< T >::ValueType