![]() |
Matrices are just as easy and intuitive to create as vectors. Still, there are a few rules to be aware of:
StaticMatrix
or HybridMatrix
are default initialized (i.e. built-in data types are initialized to 0, class types are initialized via the default constructor).DynamicMatrix
or CompressedMatrix
remain uninitialized if they are of built-in type and are default constructed if they are of class type.
The DynamicMatrix
, HybridMatrix
, and CompressedMatrix
classes offer a constructor that allows to immediately give the matrices a specific number of rows and columns:
Note that dense matrices (in this case DynamicMatrix
and HybridMatrix
) immediately allocate enough capacity for all matrix elements. Sparse matrices on the other hand (in this example CompressedMatrix
) merely acquire the size, but don't necessarily allocate memory.
All dense matrix classes offer a constructor for a direct, homogeneous initialization of all matrix elements. In contrast, for sparse matrices the predicted number of non-zero elements can be specified.
Alternatively, all dense matrix classes offer a constructor for an initialization with a dynamic or static array. If the matrix is initialized from a dynamic array, the constructor expects the dimensions of values provided by the array as first and second argument, the array as third argument. In case of a static array, the fixed size of the array is used:
In addition, all dense matrix classes can be directly initialized by means of an initializer list:
All dense and sparse matrices can be created as a copy of another dense or sparse matrix.
Note that it is not possible to create a StaticMatrix
as a copy of a matrix with a different number of rows and/or columns:
There are several types of assignment to dense and sparse matrices: Homogeneous Assignment, Array Assignment, Copy Assignment, and Compound Assignment.
It is possible to assign the same value to all elements of a dense matrix. All dense matrix classes provide an according assignment operator:
Dense matrices can also be assigned a static array:
Note that the dimensions of the static array have to match the size of a StaticMatrix
, whereas a DynamicMatrix
is resized according to the array dimensions:
Alternatively, it is possible to directly assign an initializer list to a dense matrix:
All kinds of matrices can be assigned to each other. The only restriction is that since a StaticMatrix
cannot change its size, the assigned matrix must match both in the number of rows and in the number of columns.
Compound assignment is also available for matrices: addition assignment, subtraction assignment, and multiplication assignment. In contrast to plain assignment, however, the number of rows and columns of the two operands have to match according to the arithmetic operation.
Note that the multiplication assignment potentially changes the number of columns of the target matrix:
Since a StaticMatrix
cannot change its size, only a square StaticMatrix can be used in a multiplication assignment with other square matrices of the same dimensions.
The easiest way to access a specific dense or sparse matrix element is via the function call operator. The indices to access a matrix are zero-based:
Since dense matrices allocate enough memory for all contained elements, using the function call operator on a dense matrix directly returns a reference to the accessed value. In case of a sparse matrix, if the accessed value is currently not contained in the matrix, the value is inserted into the matrix prior to returning a reference to the value, which can be much more expensive than the direct access to a dense matrix. Consider the following example:
Although the compressed matrix is only used for read access within the for loop, using the function call operator temporarily inserts 16 non-zero elements into the matrix. Therefore, all matrices (sparse as well as dense) offer an alternate way via the begin()
, cbegin()
, end()
and cend()
functions to traverse all contained elements by iterator. Note that it is not possible to traverse all elements of the matrix, but that it is only possible to traverse elements in a row/column-wise fashion. In case of a non-const matrix, begin()
and end()
return an Iterator
, which allows a manipulation of the non-zero value, in case of a constant matrix or in case cbegin()
or cend()
are used a ConstIterator
is returned:
Note that begin()
, cbegin()
, end()
, and cend()
are also available as free functions:
Whereas a dense matrix always provides enough capacity to store all matrix elements, a sparse matrix only stores the non-zero elements. Therefore it is necessary to explicitly add elements to the matrix. The first possibility to add elements to a sparse matrix is the function call operator:
In case the element at the given position is not yet contained in the sparse matrix, it is automatically inserted. Otherwise the old value is replaced by the new value 2. The operator returns a reference to the sparse vector element.
An alternative is the set()
function: In case the element is not yet contained in the matrix the element is inserted, else the element's value is modified:
However, insertion of elements can be better controlled via the insert()
function. In contrast to the function call operator and the set()
function it emits an exception in case the element is already contained in the matrix. In order to check for this case, the find()
function can be used:
Although the insert()
function is very flexible, due to performance reasons it is not suited for the setup of large sparse matrices. A very efficient, yet also very low-level way to fill a sparse matrix is the append()
function. It requires the sparse matrix to provide enough capacity to insert a new element in the specified row. Additionally, the index of the new element must be larger than the index of the previous element in the same row. Violating these conditions results in undefined behavior!
The most efficient way to fill a sparse matrix with elements, however, is a combination of reserve()
, append()
, and the finalize()
function:
The current number of rows of a matrix can be acquired via the rows()
member function:
Alternatively, the free functions rows()
can be used to query the current number of rows of a matrix. In contrast to the member function, the free function can also be used to query the number of rows of a matrix expression:
The current number of columns of a matrix can be acquired via the columns()
member function:
There is also a free function columns()
available, which can also be used to query the number of columns of a matrix expression:
The capacity()
member function returns the internal capacity of a dense or sparse matrix. Note that the capacity of a matrix doesn't have to be equal to the size of a matrix. In case of a dense matrix the capacity will always be greater or equal than the total number of elements of the matrix. In case of a sparse matrix, the capacity will usually be much less than the total number of elements.
There is also a free function capacity()
available to query the capacity. However, please note that this function cannot be used to query the capacity of a matrix expression:
For both dense and sparse matrices the current number of non-zero elements can be queried via the nonZeros()
member function. In case of matrices there are two flavors of the nonZeros()
function: One returns the total number of non-zero elements in the matrix, the second returns the number of non-zero elements in a specific row (in case of a row-major matrix) or column (in case of a column-major matrix). Sparse matrices directly return their number of non-zero elements, dense matrices traverse their elements and count the number of non-zero elements.
The free nonZeros()
function can also be used to query the number of non-zero elements in a matrix expression. However, the result is not the exact number of non-zero elements, but may be a rough estimation:
The dimensions of a StaticMatrix
are fixed at compile time by the second and third template parameter and a CustomMatrix
cannot be resized. In contrast, the number or rows and columns of DynamicMatrix
, HybridMatrix
, and CompressedMatrix
can be changed at runtime:
Note that resizing a matrix invalidates all existing views (see e.g. Submatrices) on the matrix:
When the internal capacity of a matrix is no longer sufficient, the allocation of a larger junk of memory is triggered. In order to avoid frequent reallocations, the reserve()
function can be used up front to set the internal capacity:
Additionally it is possible to reserve memory in a specific row (for a row-major matrix) or column (for a column-major matrix):
In order to reset all elements of a dense or sparse matrix, the reset()
function can be used. The number of rows and columns of the matrix are preserved:
Alternatively, only a single row or column of the matrix can be resetted:
In order to reset a row of a column-major matrix or a column of a row-major matrix, use a row or column view (see Rows and views_colums).
In order to return a matrix to its default state (i.e. the state of a default constructed matrix), the clear()
function can be used:
The isnan()
function provides the means to check a dense or sparse matrix for non-a-number elements:
If at least one element of the matrix is not-a-number, the function returns true
, otherwise it returns false
. Please note that this function only works for matrices with floating point elements. The attempt to use it for a matrix with a non-floating point element type results in a compile time error.
The isDefault()
function returns whether the given dense or sparse matrix is in default state:
A matrix is in default state if it appears to just have been default constructed. All resizable matrices (HybridMatrix
, DynamicMatrix
, or CompressedMatrix
) and CustomMatrix
are in default state if its size is equal to zero. A non-resizable matrix (StaticMatrix
and all submatrices) is in default state if all its elements are in default state. For instance, in case the matrix is instantiated for a built-in integral or floating point data type, the function returns true
in case all matrix elements are 0 and false
in case any matrix element is not 0.
Whether a dense or sparse matrix is a square matrix (i.e. if the number of rows is equal to the number of columns) can be checked via the isSquare()
function:
Via the isSymmetric()
function it is possible to check whether a dense or sparse matrix is symmetric:
Note that non-square matrices are never considered to be symmetric!
In order to check if all matrix elements are identical, the isUniform
function can be used:
Note that in case of a sparse matrix also the zero elements are also taken into account!
Via the isLower()
function it is possible to check whether a dense or sparse matrix is lower triangular:
Note that non-square matrices are never considered to be lower triangular!
Via the isUniLower()
function it is possible to check whether a dense or sparse matrix is lower unitriangular:
Note that non-square matrices are never considered to be lower unitriangular!
Via the isStrictlyLower()
function it is possible to check whether a dense or sparse matrix is strictly lower triangular:
Note that non-square matrices are never considered to be strictly lower triangular!
Via the isUpper()
function it is possible to check whether a dense or sparse matrix is upper triangular:
Note that non-square matrices are never considered to be upper triangular!
Via the isUniUpper()
function it is possible to check whether a dense or sparse matrix is upper unitriangular:
Note that non-square matrices are never considered to be upper unitriangular!
Via the isStrictlyUpper()
function it is possible to check whether a dense or sparse matrix is strictly upper triangular:
Note that non-square matrices are never considered to be strictly upper triangular!
The isDiagonal()
function checks if the given dense or sparse matrix is a diagonal matrix, i.e. if it has only elements on its diagonal and if the non-diagonal elements are default elements:
Note that non-square matrices are never considered to be diagonal!
The isIdentity()
function checks if the given dense or sparse matrix is an identity matrix, i.e. if all diagonal elements are 1 and all non-diagonal elements are 0:
Note that non-square matrices are never considered to be identity matrices!
The min()
and the max()
functions return the smallest and largest element of the given dense or sparse matrix, respectively:
In case the matrix currently has 0 rows or 0 columns, both functions return 0. Additionally, in case a given sparse matrix is not completely filled, the zero elements are taken into account. For example: the following compressed matrix has only 2 non-zero elements. However, the minimum of this matrix is 0:
Also note that the min()
and max()
functions can be used to compute the smallest and largest element of a matrix expression:
The abs()
function can be used to compute the absolute values of each element of a matrix. For instance, the following computation
results in the matrix
The floor()
and ceil()
functions can be used to round down/up each element of a matrix, respectively:
The conj()
function can be applied on a dense or sparse matrix to compute the complex conjugate of each element of the matrix:
Additionally, matrices can be conjugated in-place via the conjugate()
function:
The real()
function can be used on a dense or sparse matrix to extract the real part of each element of the matrix:
The imag()
function can be used on a dense or sparse matrix to extract the imaginary part of each element of the matrix:
Via the sqrt()
and invsqrt()
functions the (inverse) square root of each element of a matrix can be computed:
Note that in case of sparse matrices only the non-zero elements are taken into account!
The cbrt()
and invcbrt()
functions can be used to compute the the (inverse) cubic root of each element of a matrix:
Note that in case of sparse matrices only the non-zero elements are taken into account!
The clip()
function can be used to restrict all elements of a matrix to a specific range:
Note that in case of sparse matrices only the non-zero elements are taken into account!
The pow()
function can be used to compute the exponential value of each element of a matrix:
exp()
computes the base e exponential of each element of a matrix:
Note that in case of sparse matrices only the non-zero elements are taken into account!
The log()
and log10()
functions can be used to compute the natural and common logarithm of each element of a matrix:
The following trigonometric functions are available for both dense and sparse matrices:
Note that in case of sparse matrices only the non-zero elements are taken into account!
The following hyperbolic functions are available for both dense and sparse matrices:
The erf()
and erfc()
functions compute the (complementary) error function of each element of a matrix:
Note that in case of sparse matrices only the non-zero elements are taken into account!
Via the forEach()
function it is possible to execute custom operations on dense and sparse matrices. For instance, the following example demonstrates a custom square root computation via a lambda:
Although the computation can be parallelized it is not vectorized and thus cannot perform at peak performance. However, it is also possible to create vectorized custom operations. See Custom Operations for a detailed overview of the possibilities of custom operations.
Matrices can be transposed via the trans()
function. Row-major matrices are transposed into a column-major matrix and vice versa:
Additionally, matrices can be transposed in-place via the transpose()
function:
Note however that the transpose operation fails if ...
The conjugate transpose of a dense or sparse matrix (also called adjoint matrix, Hermitian conjugate, or transjugate) can be computed via the ctrans()
function:
Note that the ctrans()
function has the same effect as manually applying the conj()
and trans()
function in any order:
The ctranspose()
function can be used to perform an in-place conjugate transpose operation:
Note however that the conjugate transpose operation fails if ...
The determinant of a square dense matrix can be computed by means of the det()
function:
In case the given dense matrix is not a square matrix, a std::invalid_argument
exception is thrown.
det()
function can only be used for dense matrices with float
, double
, complex<float>
or complex<double>
element type. The attempt to call the function with matrices of any other element type or with a sparse matrix results in a compile time error!
Via the
it is possible to completely swap the contents of two matrices of the same type:swap()
function
The inverse of a square dense matrix can be computed via the inv()
function:
Alternatively, an in-place inversion of a dense matrix can be performed via the invert()
function:
Both the inv()
and the invert()
functions will automatically select the most suited matrix inversion algorithm depending on the size and type of the given matrix. For small matrices of up to 6x6, both functions use manually optimized kernels for maximum performance. For matrices larger than 6x6 the inversion is performed by means of the most suited matrix decomposition method: In case of a general matrix the LU decomposition is used, for symmetric matrices the LDLT decomposition is applied, for Hermitian matrices the LDLH decomposition is performed, and for triangular matrices the inverse is computed via a forward or back substitution.
In case the type of the matrix does not provide additional compile time information about its structure (symmetric, lower, upper, diagonal, ...), the information can be provided manually when calling the invert()
function:
Alternatively, via the invert()
function it is possible to explicitly specify the inversion algorithm:
Whereas the inversion by means of an LU decomposition works for every general square matrix, the inversion by LDLT only works for symmetric indefinite matrices, the inversion by LDLH is restricted to Hermitian indefinite matrices and the Cholesky decomposition (LLH) only works for Hermitian positive definite matrices. Please note that it is in the responsibility of the function caller to guarantee that the selected algorithm is suited for the given matrix. In case this precondition is violated the result can be wrong and might not represent the inverse of the given matrix!
For both the inv()
and invert()
function the matrix inversion fails if ...
In all failure cases either a compilation error is created if the failure can be predicted at compile time or a std::invalid_argument
exception is thrown.
float
, double
, complex<float>
or complex<double>
element type. The attempt to call the function with matrices of any other element type or with a sparse matrix results in a compile time error!inv()
function. Also, it is not possible to access individual elements via the function call operator on the expression object:
float
, double
, complex<float>
or complex<double>
element type. The attempt to call the function with matrices of any other element type or with a sparse matrix results in a compile time error!The LU decomposition of a dense matrix can be computed via the lu()
function:
The function works for both rowMajor
and columnMajor
matrices. Note, however, that the three matrices A
, L
and U
are required to have the same storage order. Also, please note that the way the permutation matrix P
needs to be applied differs between row-major and column-major matrices, since the algorithm uses column interchanges for row-major matrices and row interchanges for column-major matrices.
Furthermore, lu()
can be used with adaptors. For instance, the following example demonstrates the LU decomposition of a symmetric matrix into a lower and upper triangular matrix:
The Cholesky (LLH) decomposition of a dense matrix can be computed via the llh()
function:
The function works for both rowMajor
and columnMajor
matrices and the two matrices A
and L
can have any storage order.
Furthermore, llh()
can be used with adaptors. For instance, the following example demonstrates the LLH decomposition of a symmetric matrix into a lower triangular matrix:
The QR decomposition of a dense matrix can be computed via the qr()
function:
The function works for both rowMajor
and columnMajor
matrices and the three matrices A
, Q
and R
can have any storage order.
Furthermore, qr()
can be used with adaptors. For instance, the following example demonstrates the QR decomposition of a symmetric matrix into a general matrix and an upper triangular matrix:
Similar to the QR decomposition, the RQ decomposition of a dense matrix can be computed via the rq()
function:
The function works for both rowMajor
and columnMajor
matrices and the three matrices A
, R
and Q
can have any storage order.
Also the rq()
function can be used in combination with matrix adaptors. For instance, the following example demonstrates the RQ decomposition of an Hermitian matrix into a general matrix and an upper triangular matrix:
The QL decomposition of a dense matrix can be computed via the ql()
function:
The function works for both rowMajor
and columnMajor
matrices and the three matrices A
, Q
and L
can have any storage order.
Also the ql()
function can be used in combination with matrix adaptors. For instance, the following example demonstrates the QL decomposition of a symmetric matrix into a general matrix and a lower triangular matrix:
The LQ decomposition of a dense matrix can be computed via the lq()
function:
The function works for both rowMajor
and columnMajor
matrices and the three matrices A
, L
and Q
can have any storage order.
Furthermore, lq()
can be used with adaptors. For instance, the following example demonstrates the LQ decomposition of an Hermitian matrix into a lower triangular matrix and a general matrix:
Previous: Matrix Types Next: Adaptors