Dune Core Modules (2.10.0)

A dense n x m matrix. More...

#include <dune/common/fmatrix.hh>

Public Types

typedef Traits::derived_type derived_type
 type of derived matrix class
 
typedef Traits::value_type value_type
 export the type representing the field
 
typedef Traits::value_type field_type
 export the type representing the field
 
typedef Traits::value_type block_type
 export the type representing the components
 
typedef DenseIterator< DenseMatrix, row_type, row_referenceIterator
 Iterator class for sequential access.
 
typedef Iterator iterator
 typedef for stl compliant access
 
typedef Iterator RowIterator
 rename the iterators for easier access
 
typedef std::remove_reference< row_reference >::type::Iterator ColIterator
 rename the iterators for easier access
 
typedef DenseIterator< const DenseMatrix, const row_type, const_row_referenceConstIterator
 Iterator class for sequential access.
 
typedef ConstIterator const_iterator
 typedef for stl compliant access
 
typedef ConstIterator ConstRowIterator
 rename the iterators for easier access
 
typedef std::remove_reference< const_row_reference >::type::ConstIterator ConstColIterator
 rename the iterators for easier access
 

Public Member Functions

constexpr FieldMatrix ()=default
 Default constructor.
 
constexpr FieldMatrix (std::initializer_list< Dune::FieldVector< K, cols > > const &l)
 Constructor initializing the matrix from a list of vector.
 
FieldMatrixoperator= (const FieldMatrix &)=default
 copy assignment operator
 
template<typename T >
FieldMatrixoperator= (const FieldMatrix< T, ROWS, COLS > &x)
 copy assignment from FieldMatrix over a different field
 
template<typename T , int rows, int cols>
FieldMatrixoperator= (FieldMatrix< T, rows, cols > const &)=delete
 no copy assignment from FieldMatrix of different size
 
FieldMatrix< K, COLS, ROWS > transposed () const
 Return transposed of the matrix as FieldMatrix.
 
template<int l>
FieldMatrix< K, l, colsleftmultiplyany (const FieldMatrix< K, l, rows > &M) const
 Multiplies M from the left to this matrix, this matrix is not modified.
 
template<int r, int c>
FieldMatrixrightmultiply (const FieldMatrix< K, r, c > &M)
 Multiplies M from the right to this matrix.
 
template<int l>
FieldMatrix< K, rows, l > rightmultiplyany (const FieldMatrix< K, cols, l > &M) const
 Multiplies M from the right to this matrix, this matrix is not modified.
 
row_reference operator[] (size_type i)
 random access
 
size_type size () const
 size method (number of rows)
 
Iterator begin ()
 begin iterator
 
ConstIterator begin () const
 begin iterator
 
Iterator end ()
 end iterator
 
ConstIterator end () const
 end iterator
 
Iterator beforeEnd ()
 
ConstIterator beforeEnd () const
 
Iterator beforeBegin ()
 
ConstIterator beforeBegin () const
 
derived_typeoperator+= (const DenseMatrix< Other > &x)
 vector space addition
 
derived_type operator- () const
 Matrix negation.
 
derived_typeoperator-= (const DenseMatrix< Other > &x)
 vector space subtraction
 
derived_typeoperator*= (const field_type &k)
 vector space multiplication with scalar
 
derived_typeoperator/= (const field_type &k)
 vector space division by scalar
 
derived_typeaxpy (const field_type &a, const DenseMatrix< Other > &x)
 vector space axpy operation (*this += a x)
 
bool operator== (const DenseMatrix< Other > &x) const
 Binary matrix comparison.
 
bool operator!= (const DenseMatrix< Other > &x) const
 Binary matrix incomparison.
 
void mv (const X &x, Y &y) const
 y = A x
 
void mtv (const X &x, Y &y) const
 y = A^T x
 
void umv (const X &x, Y &y) const
 y += A x
 
void umtv (const X &x, Y &y) const
 y += A^T x
 
void umhv (const X &x, Y &y) const
 y += A^H x
 
void mmv (const X &x, Y &y) const
 y -= A x
 
void mmtv (const X &x, Y &y) const
 y -= A^T x
 
void mmhv (const X &x, Y &y) const
 y -= A^H x
 
void usmv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
 y += alpha A x
 
void usmtv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
 y += alpha A^T x
 
void usmhv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
 y += alpha A^H x
 
FieldTraits< value_type >::real_type frobenius_norm () const
 frobenius norm: sqrt(sum over squared values of entries)
 
FieldTraits< value_type >::real_type frobenius_norm2 () const
 square of frobenius norm, need for block recursion
 
FieldTraits< vt >::real_type infinity_norm () const
 infinity norm (row sum norm, how to generalize for blocks?)
 
FieldTraits< vt >::real_type infinity_norm () const
 infinity norm (row sum norm, how to generalize for blocks?)
 
FieldTraits< vt >::real_type infinity_norm_real () const
 simplified infinity norm (uses Manhattan norm for complex values)
 
FieldTraits< vt >::real_type infinity_norm_real () const
 simplified infinity norm (uses Manhattan norm for complex values)
 
void solve (V1 &x, const V2 &b, bool doPivoting=true) const
 Solve system A x = b. More...
 
void invert (bool doPivoting=true)
 Compute inverse. More...
 
field_type determinant (bool doPivoting=true) const
 calculates the determinant of this matrix
 
FieldMatrix< K, ROWS, COLS > & leftmultiply (const DenseMatrix< M2 > &M)
 Multiplies M from the left to this matrix.
 
FieldMatrix< K, ROWS, COLS > & rightmultiply (const DenseMatrix< M2 > &M)
 Multiplies M from the right to this matrix.
 
constexpr size_type N () const
 number of rows
 
constexpr size_type M () const
 number of columns
 
constexpr size_type rows () const
 number of rows
 
constexpr size_type cols () const
 number of columns
 
bool exists (size_type i, size_type j) const
 return true when (i,j) is in pattern
 

Static Public Attributes

static constexpr int rows = ROWS
 The number of rows.
 
static constexpr int cols = COLS
 The number of columns.
 
static constexpr int blocklevel
 The number of block levels we contain. This is the leaf, that is, 1.
 

Static Protected Member Functions

static void luDecomposition (DenseMatrix< FieldMatrix< K, ROWS, COLS > > &A, Func func, Mask &nonsingularLanes, bool throwEarly, bool doPivoting)
 do an LU-Decomposition on matrix A More...
 

Friends

template<class OtherScalar >
auto operator+ (const FieldMatrix &matrixA, const FieldMatrix< OtherScalar, ROWS, COLS > &matrixB)
 vector space addition – two-argument version
 
template<class OtherScalar >
auto operator- (const FieldMatrix &matrixA, const FieldMatrix< OtherScalar, ROWS, COLS > &matrixB)
 vector space subtraction – two-argument version
 
template<class Scalar , std::enable_if_t< IsNumber< Scalar >::value, int > = 0>
auto operator* (const FieldMatrix &matrix, Scalar scalar)
 vector space multiplication with scalar
 
template<class Scalar , std::enable_if_t< IsNumber< Scalar >::value, int > = 0>
auto operator* (Scalar scalar, const FieldMatrix &matrix)
 vector space multiplication with scalar
 
template<class Scalar , std::enable_if_t< IsNumber< Scalar >::value, int > = 0>
auto operator/ (const FieldMatrix &matrix, Scalar scalar)
 vector space division by scalar
 
template<class OtherScalar , int otherCols>
auto operator* (const FieldMatrix &matrixA, const FieldMatrix< OtherScalar, COLS, otherCols > &matrixB)
 Matrix-matrix multiplication.
 
template<class OtherMatrix , std::enable_if_t< Impl::IsStaticSizeMatrix_v< OtherMatrix > and not Impl::IsFieldMatrix_v< OtherMatrix >, int > = 0>
auto operator* (const FieldMatrix &matrixA, const OtherMatrix &matrixB)
 Matrix-matrix multiplication. More...
 
template<class OtherMatrix , std::enable_if_t< Impl::IsStaticSizeMatrix_v< OtherMatrix > and not Impl::IsFieldMatrix_v< OtherMatrix >, int > = 0>
auto operator* (const OtherMatrix &matrixA, const FieldMatrix &matrixB)
 Matrix-matrix multiplication. More...
 

Detailed Description

template<class K, int ROWS, int COLS>
class Dune::FieldMatrix< K, ROWS, COLS >

A dense n x m matrix.

Matrices represent linear maps from a vector space V to a vector space W. This class represents such a linear map by storing a two-dimensional array of numbers of a given field type K. The number of rows and columns is given at compile time.

Examples
recipe-integration.cc.

Member Function Documentation

◆ beforeBegin() [1/2]

Iterator Dune::DenseMatrix< FieldMatrix< K, ROWS, COLS > >::beforeBegin ( )
inlineinherited
Returns
an iterator that is positioned before the first entry of the vector.

◆ beforeBegin() [2/2]

ConstIterator Dune::DenseMatrix< FieldMatrix< K, ROWS, COLS > >::beforeBegin ( ) const
inlineinherited
Returns
an iterator that is positioned before the first entry of the vector.

◆ beforeEnd() [1/2]

Iterator Dune::DenseMatrix< FieldMatrix< K, ROWS, COLS > >::beforeEnd ( )
inlineinherited
Returns
an iterator that is positioned before the end iterator of the vector, i.e. at the last entry.

◆ beforeEnd() [2/2]

ConstIterator Dune::DenseMatrix< FieldMatrix< K, ROWS, COLS > >::beforeEnd ( ) const
inlineinherited
Returns
an iterator that is positioned before the end iterator of the vector. i.e. at the last element

◆ invert()

void Dune::DenseMatrix< FieldMatrix< K, ROWS, COLS > >::invert ( bool  doPivoting = true)
inherited

Compute inverse.

Exceptions
FMatrixErrorif the matrix is singular

◆ luDecomposition()

static void Dune::DenseMatrix< FieldMatrix< K, ROWS, COLS > >::luDecomposition ( DenseMatrix< FieldMatrix< K, ROWS, COLS > > &  A,
Func  func,
Mask &  nonsingularLanes,
bool  throwEarly,
bool  doPivoting 
)
staticprotectedinherited

do an LU-Decomposition on matrix A

Parameters
AThe matrix to decompose, and to store the result in.
funcFunctor used for swapping lanes and to conduct the elimination. Depending on the functor, luDecomposition() can be used for solving, for inverting, or to compute the determinant.
nonsingularLanesSimdMask of lanes that are nonsingular.
throwEarlyWhether to throw an FMatrixError immediately as soon as one lane is discovered to be singular. If false, do not throw, instead continue until finished or all lanes are singular, and exit via return in both cases.
doPivotingEnable pivoting.

There are two modes of operation:

  • Terminate as soon as one lane is discovered to be singular. Early termination is done by throwing an FMatrixError. On entry, Simd::allTrue(nonsingularLanes) and throwEarly==true should hold. After early termination, the contents of A should be considered bogus, and nonsingularLanes has the lane(s) that triggered the early termination unset. There may be more singular lanes than the one reported in nonsingularLanes, which just haven't been discovered yet; so the value of nonsingularLanes is mostly useful for diagnostics.
  • Terminate only when all lanes are discovered to be singular. Use this when you want to apply special postprocessing in singular lines (e.g. setting the determinant of singular lanes to 0 in determinant()). On entry, nonsingularLanes may have any value and throwEarly==false should hold. The function will not throw an exception if some lanes are discovered to be singular, instead it will continue running until all lanes are singular or until finished, and terminate only via normal return. On exit, nonsingularLanes contains the map of lanes that are valid in A.

◆ solve()

void Dune::DenseMatrix< FieldMatrix< K, ROWS, COLS > >::solve ( V1 &  x,
const V2 &  b,
bool  doPivoting = true 
) const
inherited

Solve system A x = b.

Exceptions
FMatrixErrorif the matrix is singular

Friends And Related Function Documentation

◆ operator* [1/2]

template<class K , int ROWS, int COLS>
template<class OtherMatrix , std::enable_if_t< Impl::IsStaticSizeMatrix_v< OtherMatrix > and not Impl::IsFieldMatrix_v< OtherMatrix >, int > = 0>
auto operator* ( const FieldMatrix< K, ROWS, COLS > &  matrixA,
const OtherMatrix &  matrixB 
)
friend

Matrix-matrix multiplication.

This implements multiplication of a FieldMatrix with another matrix of type OtherMatrix. The latter has to provide OtherMatrix::field_type, OtherMatrix::cols, and OtherMatrix::mtv(x,y).

◆ operator* [2/2]

template<class K , int ROWS, int COLS>
template<class OtherMatrix , std::enable_if_t< Impl::IsStaticSizeMatrix_v< OtherMatrix > and not Impl::IsFieldMatrix_v< OtherMatrix >, int > = 0>
auto operator* ( const OtherMatrix &  matrixA,
const FieldMatrix< K, ROWS, COLS > &  matrixB 
)
friend

Matrix-matrix multiplication.

This implements multiplication of another matrix of type OtherMatrix with a FieldMatrix. The former has to provide OtherMatrix::field_type, OtherMatrix::rows, and OtherMatrix::mv(x,y).


The documentation for this class was generated from the following files:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 26, 23:30, 2024)