Dune Core Modules (2.3.1)

Dune::FieldMatrix< K, ROWS, COLS > Class Template Reference

A dense n x m matrix. More...

#include <dune/common/fmatrix.hh>

Public Types

enum  { rows = ROWS , cols = COLS }
 export size More...
 
enum  
 We are at the leaf of the block recursion.
 
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_typeIterator
 Iterator class for sequential access.
 
typedef Iterator iterator
 typedef for stl compliant access
 
typedef Iterator RowIterator
 rename the iterators for easier access
 
typedef row_type::Iterator ColIterator
 rename the iterators for easier access
 
typedef DenseIterator< const DenseMatrix, const row_typeConstIterator
 Iterator class for sequential access.
 
typedef ConstIterator const_iterator
 typedef for stl compliant access
 
typedef ConstIterator ConstRowIterator
 rename the iterators for easier access
 
typedef row_type::ConstIterator ConstColIterator
 rename the iterators for easier access
 

Public Member Functions

 FieldMatrix ()
 Default constructor.
 
template<class Other >
 FieldMatrix (const Other &other)
 Constructor initializing the whole matrix with a scalar.
 
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.
 
FieldMatrixrightmultiply (const FieldMatrix< K, cols, cols > &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
 
DenseMatrixoperator+= (const DenseMatrix< Other > &y)
 vector space addition
 
DenseMatrixoperator-= (const DenseMatrix< Other > &y)
 vector space subtraction
 
DenseMatrixoperator*= (const field_type &k)
 vector space multiplication with scalar
 
DenseMatrixoperator/= (const field_type &k)
 vector space division by scalar
 
DenseMatrixaxpy (const field_type &k, const DenseMatrix< Other > &y)
 vector space axpy operation (*this += k y)
 
bool operator== (const DenseMatrix< Other > &y) const
 Binary matrix comparison.
 
bool operator!= (const DenseMatrix< Other > &y) 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 field_type &alpha, const X &x, Y &y) const
 y += alpha A x
 
void usmtv (const field_type &alpha, const X &x, Y &y) const
 y += alpha A^T x
 
void usmhv (const 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< value_type >::real_type infinity_norm () const
 infinity norm (row sum norm, how to generalize for blocks?)
 
FieldTraits< value_type >::real_type infinity_norm_real () const
 simplified infinity norm (uses Manhattan norm for complex values)
 
void solve (V &x, const V &b) const
 Solve system A x = b. More...
 
void invert ()
 Compute inverse. More...
 
field_type determinant () 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.
 
size_type N () const
 number of rows
 
size_type M () const
 number of columns
 
size_type rows () const
 number of rows
 
size_type cols () const
 number of columns
 
bool exists (size_type i, size_type j) const
 return true when (i,j) is in pattern
 

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.

Member Enumeration Documentation

◆ anonymous enum

template<class K , int ROWS, int COLS>
anonymous enum

export size

Enumerator
rows 

The number of rows.

cols 

The number of columns.

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 ( )
inherited

Compute inverse.

Exceptions
FMatrixErrorif the matrix is singular

◆ solve()

void Dune::DenseMatrix< FieldMatrix< K, ROWS, COLS > >::solve ( V &  x,
const V &  b 
) 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  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)