Dune Core Modules (2.8.0)
•All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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_type, row_reference > | Iterator |
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_reference > | ConstIterator |
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. | |
FieldMatrix (std::initializer_list< Dune::FieldVector< K, cols > > const &l) | |
Constructor initializing the matrix from a list of vector. | |
FieldMatrix & | operator= (const FieldMatrix &)=default |
copy assignment operator | |
template<typename T > | |
FieldMatrix & | operator= (const FieldMatrix< T, ROWS, COLS > &x) |
copy assignment from FieldMatrix over a different field | |
template<typename T , int rows, int cols> | |
FieldMatrix & | operator= (FieldMatrix< T, rows, cols > const &)=delete |
no copy assignment from FieldMatrix of different size | |
template<int l> | |
FieldMatrix< K, l, cols > | leftmultiplyany (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> | |
FieldMatrix & | rightmultiply (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_type & | operator+= (const DenseMatrix< Other > &x) |
vector space addition | |
derived_type | operator- () const |
Matrix negation. | |
derived_type & | operator-= (const DenseMatrix< Other > &x) |
vector space subtraction | |
derived_type & | operator*= (const field_type &k) |
vector space multiplication with scalar | |
derived_type & | operator/= (const field_type &k) |
vector space division by scalar | |
derived_type & | axpy (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 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... | |
Detailed Description
template<class K, int ROWS, int COLS>
class Dune::FieldMatrix< K, ROWS, COLS >
class Dune::FieldMatrix< K, ROWS, COLS >
A dense n x m matrix.
work around a problem of FieldMatrix/FieldVector, there is no unique way to obtain the size of a class
- Deprecated:
- VectorSize is deprecated; please call the 'size()' method directly instead. This will be removed after Dune 2.8.
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 Enumeration Documentation
◆ anonymous enum
template<class K , int ROWS, int COLS>
anonymous enum |
Member Function Documentation
◆ beforeBegin() [1/2]
|
inlineinherited |
- Returns
- an iterator that is positioned before the first entry of the vector.
◆ beforeBegin() [2/2]
|
inlineinherited |
- Returns
- an iterator that is positioned before the first entry of the vector.
◆ beforeEnd() [1/2]
|
inlineinherited |
- Returns
- an iterator that is positioned before the end iterator of the vector, i.e. at the last entry.
◆ beforeEnd() [2/2]
|
inlineinherited |
- Returns
- an iterator that is positioned before the end iterator of the vector. i.e. at the last element
◆ invert()
|
inherited |
Compute inverse.
- Exceptions
-
FMatrixError if the matrix is singular
◆ luDecomposition()
|
staticprotectedinherited |
do an LU-Decomposition on matrix A
- Parameters
-
A The matrix to decompose, and to store the result in. func Functor 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.nonsingularLanes SimdMask of lanes that are nonsingular. throwEarly Whether to throw an FMatrixError
immediately as soon as one lane is discovered to be singular. Iffalse
, do not throw, instead continue until finished or all lanes are singular, and exit via return in both cases.doPivoting Enable 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)
andthrowEarly==true
should hold. After early termination, the contents ofA
should be considered bogus, andnonsingularLanes
has the lane(s) that triggered the early termination unset. There may be more singular lanes than the one reported innonsingularLanes
, which just havent been discovered yet; so the value ofnonsingularLanes
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 andthrowEarly==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 inA
.
◆ solve()
|
inherited |
Solve system A x = b.
- Exceptions
-
FMatrixError if the matrix is singular
The documentation for this class was generated from the following files:
- dune/common/densematrix.hh
- dune/common/fmatrix.hh
