DUNE-FEM (unstable)

A dense n x m matrix. More...

#include <dune/common/fmatrix.hh>

Public Types

using size_type = typename Base::size_type
 The type used for the index access and size operation.
 
using value_type = typename Base::value_type
 The type of the elements stored in the matrix.
 
using reference = value_type &
 The type used for references to the matrix entries.
 
using const_reference = const value_type &
 The type used for const references to the matrix entries.
 
using row_type = typename Base::row_type
 The type the rows of the matrix are represented by.
 
using row_reference = typename Base::row_reference
 The type used for references to the rows of the matrix.
 
using const_row_reference = typename Base::const_row_reference
 The type used for const references to the rows of the matrix.
 
typedef Traits::derived_type derived_type
 type of derived matrix class
 
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

template<class OtherK , int otherRows, int otherCols>
requires (otherRows != ROWS || otherCols != COLS)
constexpr FieldMatrixoperator= (const FieldMatrix< OtherK, otherRows, otherCols > &)=delete
 Delete assignment from FieldMatrix of different shape.
 
constexpr FieldMatrix< K, COLS, ROWS > transposed () const
 Return transposed of the matrix as FieldMatrix.
 
constexpr row_reference operator[] (size_type i)
 random access
 
constexpr size_type size () const
 size method (number of rows)
 
constexpr Iterator begin ()
 begin iterator
 
constexpr ConstIterator begin () const
 begin iterator
 
constexpr Iterator end ()
 end iterator
 
constexpr ConstIterator end () const
 end iterator
 
constexpr Iterator beforeEnd ()
 
constexpr ConstIterator beforeEnd () const
 
constexpr Iterator beforeBegin ()
 
constexpr ConstIterator beforeBegin () const
 
constexpr derived_typeoperator+= (const DenseMatrix< Other > &x)
 vector space addition
 
constexpr derived_type operator- () const
 Matrix negation.
 
constexpr derived_typeoperator-= (const DenseMatrix< Other > &x)
 vector space subtraction
 
constexpr derived_typeoperator*= (const field_type &k)
 vector space multiplication with scalar
 
constexpr derived_typeoperator/= (const field_type &k)
 vector space division by scalar
 
constexpr derived_typeaxpy (const field_type &a, const DenseMatrix< Other > &x)
 vector space axpy operation (*this += a x)
 
constexpr bool operator== (const DenseMatrix< Other > &x) const
 Binary matrix comparison.
 
constexpr bool operator!= (const DenseMatrix< Other > &x) const
 Binary matrix incomparison.
 
constexpr void mv (const X &x, Y &y) const
 y = A x
 
constexpr void mtv (const X &x, Y &y) const
 y = A^T x
 
constexpr void umv (const X &x, Y &y) const
 y += A x
 
constexpr void umtv (const X &x, Y &y) const
 y += A^T x
 
constexpr void umhv (const X &x, Y &y) const
 y += A^H x
 
constexpr void mmv (const X &x, Y &y) const
 y -= A x
 
constexpr void mmtv (const X &x, Y &y) const
 y -= A^T x
 
constexpr void mmhv (const X &x, Y &y) const
 y -= A^H x
 
constexpr void usmv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
 y += alpha A x
 
constexpr void usmtv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
 y += alpha A^T x
 
constexpr void usmhv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const
 y += alpha A^H x
 
constexpr FieldTraits< value_type >::real_type frobenius_norm () const
 frobenius norm: sqrt(sum over squared values of entries)
 
constexpr FieldTraits< value_type >::real_type frobenius_norm2 () const
 square of frobenius norm, need for block recursion
 
constexpr FieldTraits< vt >::real_type infinity_norm () const
 infinity norm (row sum norm, how to generalize for blocks?)
 
constexpr FieldTraits< vt >::real_type infinity_norm () const
 infinity norm (row sum norm, how to generalize for blocks?)
 
constexpr FieldTraits< vt >::real_type infinity_norm_real () const
 simplified infinity norm (uses Manhattan norm for complex values)
 
constexpr 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
 
constexpr bool exists (size_type i, size_type j) const
 return true when (i,j) is in pattern
 
Constructors
constexpr FieldMatrix () noexcept(std::is_nothrow_default_constructible_v< K >)
 Default constructor, making value-initialized matrix with all components set to zero.
 
constexpr FieldMatrix (std::initializer_list< Dune::FieldVector< K, cols > > const &l)
 Constructor initializing the matrix from a nested list of values.
 
template<class OtherK , int otherRows, int otherCols>
requires (otherRows != ROWS || otherCols != COLS)
constexpr FieldMatrix (const FieldMatrix< OtherK, otherRows, otherCols > &)=delete
 Delete assignment from FieldMatrix of different shape.
 
template<class OtherMatrix >
requires (HasDenseMatrixAssigner<FieldMatrix, OtherMatrix>::value)
constexpr FieldMatrix (const OtherMatrix &rhs)
 copy constructor from assignable type OtherMatrix
 
Element access
constexpr row_reference mat_access (size_type i)
 
constexpr const_row_reference mat_access (size_type i) const
 
constexpr operator const_reference () const noexcept
 Conversion operator.
 
constexpr operator reference () noexcept
 Conversion operator.
 

Static Public Member Functions

Capacity
static constexpr size_type mat_rows ()
 
static constexpr size_type mat_cols ()
 

Static Public Attributes

static constexpr std::integral_constant< int, ROWS > rows = {}
 The number of rows.
 
static constexpr std::integral_constant< 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...
 

Vector space operations

template<Concept::Number S>
requires (ROWS*COLS == 1)
constexpr FieldMatrixoperator+= (const S &scalar)
 add scalar
 
template<Concept::Number S>
requires (ROWS*COLS == 1)
constexpr FieldMatrixoperator-= (const S &scalar)
 subtract scalar
 
template<Concept::Number S>
requires (ROWS*COLS == 1)
constexpr FieldMatrixoperator*= (const S &scalar)
 multiplication with scalar
 
template<Concept::Number S>
requires (ROWS*COLS == 1)
constexpr FieldMatrixoperator/= (const S &scalar)
 division by scalar
 

Matrix-matrix multiplication

template<int l>
constexpr 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>
constexpr FieldMatrixrightmultiply (const FieldMatrix< K, r, c > &M)
 Multiplies M from the right to this matrix.
 
template<int l>
constexpr 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.
 

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]

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

◆ beforeBegin() [2/2]

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

◆ beforeEnd() [1/2]

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

◆ beforeEnd() [2/2]

constexpr ConstIterator Dune::DenseMatrix< FieldMatrix< K, ROWS, COLS > >::beforeEnd ( ) const
inlineconstexprinherited
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

The documentation for this class was generated from the following files:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Feb 16, 23:40, 2026)