dune-common
2.8.0
|
A dense n x m matrix. More...
#include <dune/common/fmatrix.hh>
Public Types | |
enum | { rows = ROWS , cols = COLS } |
export size More... | |
typedef Base::size_type | size_type |
typedef Base::row_type | row_type |
typedef Base::row_reference | row_reference |
typedef Base::const_row_reference | const_row_reference |
enum | |
We are at the leaf of the block recursion. More... | |
typedef Traits::derived_type | derived_type |
type of derived matrix class More... | |
typedef Traits::value_type | value_type |
export the type representing the field More... | |
typedef Traits::value_type | field_type |
export the type representing the field More... | |
typedef Traits::value_type | block_type |
export the type representing the components More... | |
typedef DenseIterator< DenseMatrix, row_type, row_reference > | Iterator |
Iterator class for sequential access. More... | |
typedef Iterator | iterator |
typedef for stl compliant access More... | |
typedef Iterator | RowIterator |
rename the iterators for easier access More... | |
typedef std::remove_reference< row_reference >::type::Iterator | ColIterator |
rename the iterators for easier access More... | |
typedef DenseIterator< const DenseMatrix, const row_type, const_row_reference > | ConstIterator |
Iterator class for sequential access. More... | |
typedef ConstIterator | const_iterator |
typedef for stl compliant access More... | |
typedef ConstIterator | ConstRowIterator |
rename the iterators for easier access More... | |
typedef std::remove_reference< const_row_reference >::type::ConstIterator | ConstColIterator |
rename the iterators for easier access More... | |
Public Member Functions | |
constexpr | FieldMatrix ()=default |
Default constructor. More... | |
FieldMatrix (std::initializer_list< Dune::FieldVector< K, cols > > const &l) | |
Constructor initializing the matrix from a list of vector. More... | |
template<class T , typename = std::enable_if_t<HasDenseMatrixAssigner<FieldMatrix, T>::value>> | |
FieldMatrix (T const &rhs) | |
FieldMatrix & | operator= (const FieldMatrix &)=default |
copy assignment operator More... | |
template<typename T > | |
FieldMatrix & | operator= (const FieldMatrix< T, ROWS, COLS > &x) |
copy assignment from FieldMatrix over a different field More... | |
template<typename T , int rows, int cols> | |
FieldMatrix & | operator= (FieldMatrix< T, rows, cols > const &)=delete |
no copy assignment from FieldMatrix of different size More... | |
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. More... | |
template<int r, int c> | |
FieldMatrix & | rightmultiply (const FieldMatrix< K, r, c > &M) |
Multiplies M from the right to this matrix. More... | |
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. More... | |
row_reference | mat_access (size_type i) |
const_row_reference | mat_access (size_type i) const |
row_reference | operator[] (size_type i) |
random access More... | |
const_row_reference | operator[] (size_type i) const |
size_type | size () const |
size method (number of rows) More... | |
Iterator | begin () |
begin iterator More... | |
ConstIterator | begin () const |
begin iterator More... | |
Iterator | end () |
end iterator More... | |
ConstIterator | end () const |
end iterator More... | |
Iterator | beforeEnd () |
ConstIterator | beforeEnd () const |
Iterator | beforeBegin () |
ConstIterator | beforeBegin () const |
derived_type & | operator+= (const DenseMatrix< Other > &x) |
vector space addition More... | |
derived_type | operator- () const |
Matrix negation. More... | |
derived_type & | operator-= (const DenseMatrix< Other > &x) |
vector space subtraction More... | |
derived_type & | operator*= (const field_type &k) |
vector space multiplication with scalar More... | |
derived_type & | operator/= (const field_type &k) |
vector space division by scalar More... | |
derived_type & | axpy (const field_type &a, const DenseMatrix< Other > &x) |
vector space axpy operation (*this += a x) More... | |
bool | operator== (const DenseMatrix< Other > &x) const |
Binary matrix comparison. More... | |
bool | operator!= (const DenseMatrix< Other > &x) const |
Binary matrix incomparison. More... | |
void | mv (const X &x, Y &y) const |
y = A x More... | |
void | mtv (const X &x, Y &y) const |
y = A^T x More... | |
void | umv (const X &x, Y &y) const |
y += A x More... | |
void | umtv (const X &x, Y &y) const |
y += A^T x More... | |
void | umhv (const X &x, Y &y) const |
y += A^H x More... | |
void | mmv (const X &x, Y &y) const |
y -= A x More... | |
void | mmtv (const X &x, Y &y) const |
y -= A^T x More... | |
void | mmhv (const X &x, Y &y) const |
y -= A^H x More... | |
void | usmv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const |
y += alpha A x More... | |
void | usmtv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const |
y += alpha A^T x More... | |
void | usmhv (const typename FieldTraits< Y >::field_type &alpha, const X &x, Y &y) const |
y += alpha A^H x More... | |
FieldTraits< value_type >::real_type | frobenius_norm () const |
frobenius norm: sqrt(sum over squared values of entries) More... | |
FieldTraits< value_type >::real_type | frobenius_norm2 () const |
square of frobenius norm, need for block recursion More... | |
FieldTraits< vt >::real_type | infinity_norm () const |
infinity norm (row sum norm, how to generalize for blocks?) More... | |
FieldTraits< vt >::real_type | infinity_norm () const |
infinity norm (row sum norm, how to generalize for blocks?) More... | |
FieldTraits< vt >::real_type | infinity_norm_real () const |
simplified infinity norm (uses Manhattan norm for complex values) More... | |
FieldTraits< vt >::real_type | infinity_norm_real () const |
simplified infinity norm (uses Manhattan norm for complex values) More... | |
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 More... | |
FieldMatrix< K, ROWS, COLS > & | leftmultiply (const DenseMatrix< M2 > &M) |
Multiplies M from the left to this matrix. More... | |
FieldMatrix< K, ROWS, COLS > & | rightmultiply (const DenseMatrix< M2 > &M) |
Multiplies M from the right to this matrix. More... | |
constexpr size_type | N () const |
number of rows More... | |
constexpr size_type | M () const |
number of columns More... | |
constexpr size_type | rows () const |
number of rows More... | |
constexpr size_type | cols () const |
number of columns More... | |
bool | exists ([[maybe_unused]] size_type i,[[maybe_unused]] size_type j) const |
return true when (i,j) is in pattern More... | |
Static Public Member Functions | |
static constexpr size_type | mat_rows () |
static constexpr size_type | mat_cols () |
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... | |
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
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.
|
inherited |
export the type representing the components
|
inherited |
rename the iterators for easier access
|
inherited |
typedef for stl compliant access
typedef Base::const_row_reference Dune::FieldMatrix< K, ROWS, COLS >::const_row_reference |
|
inherited |
rename the iterators for easier access
|
inherited |
Iterator class for sequential access.
|
inherited |
rename the iterators for easier access
|
inherited |
type of derived matrix class
|
inherited |
export the type representing the field
|
inherited |
Iterator class for sequential access.
|
inherited |
typedef for stl compliant access
typedef Base::row_reference Dune::FieldMatrix< K, ROWS, COLS >::row_reference |
typedef Base::row_type Dune::FieldMatrix< K, ROWS, COLS >::row_type |
|
inherited |
rename the iterators for easier access
typedef Base::size_type Dune::FieldMatrix< K, ROWS, COLS >::size_type |
|
inherited |
export the type representing the field
anonymous enum |
|
inherited |
We are at the leaf of the block recursion.
|
constexprdefault |
Default constructor.
|
inline |
Constructor initializing the matrix from a list of vector.
|
inline |
|
inlineinherited |
vector space axpy operation (*this += a x)
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
begin iterator
|
inlineinherited |
begin iterator
|
inlineconstexprinherited |
number of columns
|
inherited |
calculates the determinant of this matrix
|
inlineinherited |
end iterator
|
inlineinherited |
end iterator
|
inlineinherited |
return true when (i,j) is in pattern
|
inlineinherited |
frobenius norm: sqrt(sum over squared values of entries)
|
inlineinherited |
square of frobenius norm, need for block recursion
|
inlineinherited |
infinity norm (row sum norm, how to generalize for blocks?)
|
inlineinherited |
infinity norm (row sum norm, how to generalize for blocks?)
|
inlineinherited |
simplified infinity norm (uses Manhattan norm for complex values)
|
inlineinherited |
simplified infinity norm (uses Manhattan norm for complex values)
|
inherited |
Compute inverse.
FMatrixError | if the matrix is singular |
|
inlineinherited |
Multiplies M from the left to this matrix.
|
inline |
Multiplies M from the left to this matrix, this matrix is not modified.
|
staticprotectedinherited |
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 havent 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
.
|
inlineconstexprinherited |
number of columns
|
inline |
|
inline |
|
inlinestaticconstexpr |
|
inlinestaticconstexpr |
|
inlineinherited |
y -= A^H x
|
inlineinherited |
y -= A^T x
|
inlineinherited |
y -= A x
|
inlineinherited |
y = A^T x
|
inlineinherited |
y = A x
|
inlineconstexprinherited |
number of rows
|
inlineinherited |
Binary matrix incomparison.
|
inlineinherited |
vector space multiplication with scalar
|
inlineinherited |
vector space addition
|
inlineinherited |
Matrix negation.
|
inlineinherited |
vector space subtraction
|
inlineinherited |
vector space division by scalar
|
default |
copy assignment operator
|
inline |
copy assignment from FieldMatrix over a different field
|
delete |
no copy assignment from FieldMatrix of different size
|
inlineinherited |
Binary matrix comparison.
|
inlineinherited |
random access
|
inlineinherited |
|
inlineinherited |
Multiplies M from the right to this matrix.
|
inline |
Multiplies M from the right to this matrix.
|
inline |
Multiplies M from the right to this matrix, this matrix is not modified.
|
inlineconstexprinherited |
number of rows
|
inlineinherited |
size method (number of rows)
|
inherited |
Solve system A x = b.
FMatrixError | if the matrix is singular |
|
inlineinherited |
y += A^H x
|
inlineinherited |
y += A^T x
|
inlineinherited |
y += A x
|
inlineinherited |
y += alpha A^H x
|
inlineinherited |
y += alpha A^T x
|
inlineinherited |
y += alpha A x