dune-common 2.9.0
|
A dense n x m matrix. More...
#include <dune/common/densematrix.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 Traits::size_type | size_type |
The type used for the index access and size operation. | |
typedef Traits::row_type | row_type |
The type used to represent a row (must fulfill the Dune::DenseVector interface) | |
typedef Traits::row_reference | row_reference |
The type used to represent a reference to a row (usually row_type &) | |
typedef Traits::const_row_reference | const_row_reference |
The type used to represent a reference to a constant row (usually const row_type &) | |
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 | |
row_reference | operator[] (size_type i) |
random access | |
const_row_reference | operator[] (size_type i) const |
size_type | size () const |
size method (number of rows) | |
Iterator | begin () |
begin iterator | |
Iterator | end () |
end iterator | |
Iterator | beforeEnd () |
Iterator | beforeBegin () |
ConstIterator | begin () const |
begin iterator | |
ConstIterator | end () const |
end iterator | |
ConstIterator | beforeEnd () const |
ConstIterator | beforeBegin () const |
template<class RHS , class = std::enable_if_t< HasDenseMatrixAssigner< MAT, RHS >::value >> | |
derived_type & | operator= (const RHS &rhs) |
template<class Other > | |
derived_type & | operator+= (const DenseMatrix< Other > &x) |
vector space addition | |
derived_type | operator- () const |
Matrix negation. | |
template<class Other > | |
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 | |
template<class Other > | |
derived_type & | axpy (const field_type &a, const DenseMatrix< Other > &x) |
vector space axpy operation (*this += a x) | |
template<class Other > | |
bool | operator== (const DenseMatrix< Other > &x) const |
Binary matrix comparison. | |
template<class Other > | |
bool | operator!= (const DenseMatrix< Other > &x) const |
Binary matrix incomparison. | |
template<class X , class Y > | |
void | mv (const X &x, Y &y) const |
y = A x | |
template<class X , class Y > | |
void | mtv (const X &x, Y &y) const |
y = A^T x | |
template<class X , class Y > | |
void | umv (const X &x, Y &y) const |
y += A x | |
template<class X , class Y > | |
void | umtv (const X &x, Y &y) const |
y += A^T x | |
template<class X , class Y > | |
void | umhv (const X &x, Y &y) const |
y += A^H x | |
template<class X , class Y > | |
void | mmv (const X &x, Y &y) const |
y -= A x | |
template<class X , class Y > | |
void | mmtv (const X &x, Y &y) const |
y -= A^T x | |
template<class X , class Y > | |
void | mmhv (const X &x, Y &y) const |
y -= A^H x | |
template<class X , class Y > | |
void | usmv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const |
y += alpha A x | |
template<class X , class Y > | |
void | usmtv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const |
y += alpha A^T x | |
template<class X , class Y > | |
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 | |
template<typename vt = value_type, typename std::enable_if<!HasNaN< vt >::value, int >::type = 0> | |
FieldTraits< vt >::real_type | infinity_norm () const |
infinity norm (row sum norm, how to generalize for blocks?) | |
template<typename vt = value_type, typename std::enable_if<!HasNaN< vt >::value, int >::type = 0> | |
FieldTraits< vt >::real_type | infinity_norm_real () const |
simplified infinity norm (uses Manhattan norm for complex values) | |
template<typename vt = value_type, typename std::enable_if< HasNaN< vt >::value, int >::type = 0> | |
FieldTraits< vt >::real_type | infinity_norm () const |
infinity norm (row sum norm, how to generalize for blocks?) | |
template<typename vt = value_type, typename std::enable_if< HasNaN< vt >::value, int >::type = 0> | |
FieldTraits< vt >::real_type | infinity_norm_real () const |
simplified infinity norm (uses Manhattan norm for complex values) | |
template<class V1 , class V2 > | |
void | solve (V1 &x, const V2 &b, bool doPivoting=true) const |
Solve system A x = b. | |
void | invert (bool doPivoting=true) |
Compute inverse. | |
field_type | determinant (bool doPivoting=true) const |
calculates the determinant of this matrix | |
template<typename M2 > | |
MAT & | leftmultiply (const DenseMatrix< M2 > &M) |
Multiplies M from the left to this matrix. | |
template<typename M2 > | |
MAT & | 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 | blocklevel = 1 |
The number of block levels we contain. This is the leaf, that is, 1. | |
Static Protected Member Functions | |
template<class Func , class Mask > | |
static void | luDecomposition (DenseMatrix< MAT > &A, Func func, Mask &nonsingularLanes, bool throwEarly, bool doPivoting) |
do an LU-Decomposition on matrix A | |
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.
MAT | type of the matrix implementation |
typedef Traits::value_type Dune::DenseMatrix< MAT >::block_type |
export the type representing the components
typedef std::remove_reference<row_reference>::type::Iterator Dune::DenseMatrix< MAT >::ColIterator |
rename the iterators for easier access
typedef ConstIterator Dune::DenseMatrix< MAT >::const_iterator |
typedef for stl compliant access
typedef Traits::const_row_reference Dune::DenseMatrix< MAT >::const_row_reference |
The type used to represent a reference to a constant row (usually const row_type &)
typedef std::remove_reference<const_row_reference>::type::ConstIterator Dune::DenseMatrix< MAT >::ConstColIterator |
rename the iterators for easier access
typedef DenseIterator<const DenseMatrix,const row_type,const_row_reference> Dune::DenseMatrix< MAT >::ConstIterator |
Iterator class for sequential access.
typedef ConstIterator Dune::DenseMatrix< MAT >::ConstRowIterator |
rename the iterators for easier access
typedef Traits::derived_type Dune::DenseMatrix< MAT >::derived_type |
type of derived matrix class
typedef Traits::value_type Dune::DenseMatrix< MAT >::field_type |
export the type representing the field
typedef DenseIterator<DenseMatrix,row_type,row_reference> Dune::DenseMatrix< MAT >::Iterator |
Iterator class for sequential access.
typedef Iterator Dune::DenseMatrix< MAT >::iterator |
typedef for stl compliant access
typedef Traits::row_reference Dune::DenseMatrix< MAT >::row_reference |
The type used to represent a reference to a row (usually row_type &)
typedef Traits::row_type Dune::DenseMatrix< MAT >::row_type |
The type used to represent a row (must fulfill the Dune::DenseVector interface)
typedef Iterator Dune::DenseMatrix< MAT >::RowIterator |
rename the iterators for easier access
typedef Traits::size_type Dune::DenseMatrix< MAT >::size_type |
The type used for the index access and size operation.
typedef Traits::value_type Dune::DenseMatrix< MAT >::value_type |
export the type representing the field
|
inline |
vector space axpy operation (*this += a x)
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
begin iterator
|
inline |
begin iterator
|
inlineconstexpr |
number of columns
field_type Dune::DenseMatrix< MAT >::determinant | ( | bool | doPivoting = true | ) | const |
calculates the determinant of this matrix
|
inline |
end iterator
|
inline |
end iterator
|
inline |
return true when (i,j) is in pattern
|
inline |
frobenius norm: sqrt(sum over squared values of entries)
|
inline |
square of frobenius norm, need for block recursion
|
inline |
infinity norm (row sum norm, how to generalize for blocks?)
|
inline |
infinity norm (row sum norm, how to generalize for blocks?)
|
inline |
simplified infinity norm (uses Manhattan norm for complex values)
|
inline |
simplified infinity norm (uses Manhattan norm for complex values)
void Dune::DenseMatrix< MAT >::invert | ( | bool | doPivoting = true | ) |
Compute inverse.
FMatrixError | if the matrix is singular |
|
inline |
Multiplies M from the left to this matrix.
|
staticprotected |
do an LU-Decomposition on matrix A
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. If false , 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:
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. 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
.
|
inlineconstexpr |
number of columns
|
inline |
y -= A^H x
|
inline |
y -= A^T x
|
inline |
y -= A x
|
inline |
y = A^T x
|
inline |
y = A x
|
inlineconstexpr |
number of rows
|
inline |
Binary matrix incomparison.
|
inline |
vector space multiplication with scalar
|
inline |
vector space addition
|
inline |
Matrix negation.
|
inline |
vector space subtraction
|
inline |
vector space division by scalar
|
inline |
|
inline |
Binary matrix comparison.
|
inline |
random access
|
inline |
|
inline |
Multiplies M from the right to this matrix.
|
inlineconstexpr |
number of rows
|
inline |
size method (number of rows)
void Dune::DenseMatrix< MAT >::solve | ( | V1 & | x, |
const V2 & | b, | ||
bool | doPivoting = true |
||
) | const |
Solve system A x = b.
FMatrixError | if the matrix is singular |
|
inline |
y += A^H x
|
inline |
y += A^T x
|
inline |
y += A x
|
inline |
y += alpha A^H x
|
inline |
y += alpha A^T x
|
inline |
y += alpha A x
|
staticconstexpr |
The number of block levels we contain. This is the leaf, that is, 1.