The Eigen library is divided in a Core module and several additional modules. Each module has a corresponding header file which has to be included in order to use the module. The Dense
and Eigen
header files are provided to conveniently gain access to several modules at once.
Module | Header file | Contents |
---|---|---|
Core | #include <Eigen/Core>
| Matrix and Array classes, basic linear algebra (including triangular and selfadjoint products), array manipulation |
Geometry | #include <Eigen/Geometry>
| Transform, Translation, Scaling, Rotation2D and 3D rotations (Quaternion, AngleAxis) |
LU | Inverse, determinant, LU decompositions with solver (FullPivLU, PartialPivLU) | |
Cholesky | #include <Eigen/Cholesky>
| LLT and LDLT Cholesky factorization with solver |
Householder | #include <Eigen/Householder>
| Householder transformations; this module is used by several linear algebra modules |
SVD | SVD decompositions with least-squares solver (JacobiSVD, BDCSVD) | |
QR | QR decomposition with solver (HouseholderQR, ColPivHouseholderQR, FullPivHouseholderQR) | |
Eigenvalues | #include <Eigen/Eigenvalues>
| Eigenvalue, eigenvector decompositions (EigenSolver, SelfAdjointEigenSolver, ComplexEigenSolver) |
Sparse | #include <Eigen/Sparse>
| Sparse matrix storage and related basic linear algebra (SparseMatrix, SparseVector) (see Quick reference guide for sparse matrices for details on sparse modules) |
#include <Eigen/Dense>
| Includes Core, Geometry, LU, Cholesky, SVD, QR, and Eigenvalues header files | |
#include <Eigen/Eigen>
| Includes Dense and Sparse header files (the whole Eigen library) |
Recall: Eigen provides two kinds of dense objects: mathematical matrices and vectors which are both represented by the template class Matrix, and general 1D and 2D arrays represented by the template class Array:
Scalar
is the scalar type of the coefficients (e.g., float
, double
, bool
, int
, etc.). RowsAtCompileTime
and ColsAtCompileTime
are the number of rows and columns of the matrix as known at compile-time or Dynamic
. Options
can be ColMajor
or RowMajor
, default is ColMajor
. (see class Matrix for more options)All combinations are allowed: you can have a matrix with a fixed number of rows and a dynamic number of columns, etc. The following are all valid:
In most cases, you can simply use one of the convenience typedefs for matrices and arrays. Some examples:
Matrices | Arrays |
---|---|
Matrix<float,Dynamic,Dynamic> <=> MatrixXf
Matrix<double,Dynamic,1> <=> VectorXd
Matrix<int,1,Dynamic> <=> RowVectorXi
Matrix<float,3,3> <=> Matrix3f
Matrix<float,4,1> <=> Vector4f
Matrix< float, Dynamic, Dynamic > MatrixXf Dynamic×Dynamic matrix of type float. Definition: Matrix.h:501 |
Conversion between the matrix and array worlds:
In the rest of this document we will use the following symbols to emphasize the features which are specifics to a given kind of object:
1D objects | 2D objects | Notes | |
---|---|---|---|
Constructors | Vector4d v4;
VectorXf v5; // empty object
Map< RowVectorXf > v2(M2.data(), M2.size()) M1<< 1, 2, 3, 4, 5, 6, 7, 8, 9;Map< RowVectorXf > v1(M1.data(), M1.size()) | By default, the coefficients are left uninitialized | |
Comma initializer | |||
Comma initializer (bis) | MatrixXf::Zero(3,cols-3),
MatrixXf::Zero(rows-3,3),
cout << m;
static const IdentityReturnType Identity() Definition: CwiseNullaryOp.h:828 | output: 1 2 3 0 0 4 5 6 0 0 7 8 9 0 0 0 0 0 1 0 0 0 0 0 1 | |
Runtime info | vector.size();
vector.innerStride();
vector.data();
| matrix.rows(); matrix.cols();
matrix.innerSize(); matrix.outerSize();
matrix.innerStride(); matrix.outerStride();
matrix.data();
| Inner/Outer* are storage order dependent |
Compile-time info | ObjectType::Scalar ObjectType::RowsAtCompileTime
ObjectType::RealScalar ObjectType::ColsAtCompileTime
ObjectType::Index ObjectType::SizeAtCompileTime
| ||
Resizing | matrix.resize(nb_rows, nb_cols);
matrix.resize(Eigen::NoChange, nb_cols);
matrix.resize(nb_rows, Eigen::NoChange);
matrix.resizeLike(other_matrix);
matrix.conservativeResize(nb_rows, nb_cols);
| no-op if the new sizes match, | |
Coeff access with range checking | Range checking is disabled if | ||
Coeff access without range checking | |||
Assignment/copy | object = expression;
object_of_float = expression_of_double.cast<float>();
| the destination is automatically resized (if possible) |
Fixed-size matrix or vector | Dynamic-size matrix | Dynamic-size vector |
---|---|---|
Identity and basis vectors * | ||
x = FixedXD::Identity();
x.setIdentity();
Vector3f::UnitX() // 1 0 0
Vector3f::UnitY() // 0 1 0
Vector3f::UnitZ() // 0 0 1
| N/A
VectorXf::Unit(4,1) == Vector4f(0,1,0,0)
== Vector4f::UnitY()
static const BasisReturnType Unit(Index size, Index i) Definition: CwiseNullaryOp.h:931 |
Note that it is allowed to call any of the set*
functions to a dynamic-sized vector or matrix without passing new sizes. For instance:
Contiguous memory | float data[] = {1,2,3,4};
|
Typical usage of strides | float data[] = {1,2,3,4,5,6,7,8,9};
|
add subtract | |
scalar product | |
matrix/vector products * | |
transposition adjoint * | |
dot product inner product * | scalar = col1.adjoint() * col2;
scalar = (col1.adjoint() * col2).value();
RowVectorXd vec1(3) ScalarBinaryOpTraits< typename internal::traits< Derived >::Scalar, typename internal::traits< OtherDerived >::Scalar >::ReturnType dot(const MatrixBase< OtherDerived > &other) const Definition: Dot.h:69 |
outer product * | |
norm normalization * | |
cross product * | #include <Eigen/Geometry>
v3c = v3a.cross(v3b); // size-3 vectors
scalar = v2a.cross(v2b); // size-2 vectors
|
In addition to the aforementioned operators, Eigen supports numerous coefficient-wise operator and functions. Most of them unambiguously makes sense in array-world*. The following operators are readily available for arrays, or available through .array() for vectors and matrices:
Arithmetic operators | array1 * array2 array1 / array2 array1 *= array2 array1 /= array2
array1 + scalar array1 - scalar array1 += scalar array1 -= scalar
|
Comparisons | array1 < array2 array1 > array2 array1 < scalar array1 > scalar
array1 <= array2 array1 >= array2 array1 <= scalar array1 >= scalar
array1 == array2 array1 != array2 array1 == scalar array1 != scalar
array1.min(array2) array1.max(array2) array1.min(scalar) array1.max(scalar)
|
Trigo, power, and misc functions and the STL-like variants | array1.abs2()
array1.abs() abs(array1)
array1.sqrt() sqrt(array1)
array1.log() log(array1)
array1.log10() log10(array1)
array1.exp() exp(array1)
array1.pow(array2) pow(array1,array2)
array1.pow(scalar) pow(array1,scalar)
pow(scalar,array2)
array1.square()
array1.cube()
array1.inverse()
array1.sin() sin(array1)
array1.cos() cos(array1)
array1.tan() tan(array1)
array1.asin() asin(array1)
array1.acos() acos(array1)
array1.atan() atan(array1)
array1.sinh() sinh(array1)
array1.cosh() cosh(array1)
array1.tanh() tanh(array1)
array1.arg() arg(array1)
array1.floor() floor(array1)
array1.ceil() ceil(array1)
array1.round() round(aray1)
array1.isFinite() isfinite(array1)
array1.isInf() isinf(array1)
array1.isNaN() isnan(array1)
bfloat16 pow(const bfloat16 &a, const bfloat16 &b) Definition: BFloat16.h:626 const Eigen::CwiseUnaryOp< Eigen::internal::scalar_tanh_op< typename Derived::Scalar >, const Derived > tanh(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isinf_op< typename Derived::Scalar >, const Derived > isinf(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_atan_op< typename Derived::Scalar >, const Derived > atan(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_log10_op< typename Derived::Scalar >, const Derived > log10(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_cosh_op< typename Derived::Scalar >, const Derived > cosh(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_tan_op< typename Derived::Scalar >, const Derived > tan(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_acos_op< typename Derived::Scalar >, const Derived > acos(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_asin_op< typename Derived::Scalar >, const Derived > asin(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isnan_op< typename Derived::Scalar >, const Derived > isnan(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_arg_op< typename Derived::Scalar >, const Derived > arg(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_ceil_op< typename Derived::Scalar >, const Derived > ceil(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isfinite_op< typename Derived::Scalar >, const Derived > isfinite(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_floor_op< typename Derived::Scalar >, const Derived > floor(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_cos_op< typename Derived::Scalar >, const Derived > cos(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_round_op< typename Derived::Scalar >, const Derived > round(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_exp_op< typename Derived::Scalar >, const Derived > exp(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sin_op< typename Derived::Scalar >, const Derived > sin(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_log_op< typename Derived::Scalar >, const Derived > log(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sinh_op< typename Derived::Scalar >, const Derived > sinh(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x) |
The following coefficient-wise operators are available for all kind of expressions (matrices, vectors, and arrays), and for both real or complex scalar types:
Eigen's API | STL-like APIs* | Comments |
---|---|---|
real(array1)
imag(array1)
conj(array1)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_imag_op< typename Derived::Scalar >, const Derived > imag(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x) const Eigen::CwiseUnaryOp< Eigen::internal::scalar_real_op< typename Derived::Scalar >, const Derived > real(const Eigen::ArrayBase< Derived > &x) | // read-write, no-op for real expressions
// read-only for real, read-write for complexes
// no-op for real expressions
|
Some coefficient-wise operators are readily available for for matrices and vectors through the following cwise* methods:
Matrix API * | Via Array conversions |
---|---|
mat1.cwiseAbs2()
mat1.cwiseAbs()
mat1.cwiseSqrt()
mat1.cwiseInverse()
mat1.cwiseProduct(mat2)
mat1.cwiseQuotient(mat2)
mat1.cwiseNotEqual(mat2)
const SparseMatrixBase< OtherDerived >::template CwiseProductDenseReturnType< Derived >::Type cwiseProduct(const SparseMatrixBase< OtherDerived > &other) const Definition: MatrixBase.h:457 |
The main difference between the two API is that the one based on cwise* methods returns an expression in the matrix world, while the second one (based on .array()) returns an array expression. Recall that .array() has no cost, it only changes the available API and interpretation of the data.
It is also very simple to apply any user defined function foo
using DenseBase::unaryExpr together with std::ptr_fun (c++03, deprecated or removed in newer C++ versions), std::ref (c++11), or lambdas (c++11):
Please note that it's not possible to pass a raw function pointer to unaryExpr
, so please warp it as shown above.
Eigen provides several reduction methods such as: minCoeff() , maxCoeff() , sum() , prod() , trace() *, norm() *, squaredNorm() *, all() , and any() . All reduction operations can be done matrix-wise, column-wise or row-wise . Usage example:
internal::traits< Derived >::Scalar minCoeff() const Definition: Redux.h:518 | 1
| |
const MinCoeffReturnType minCoeff() const Definition: VectorwiseOp.h:377 | 2 3 1
| |
1
2
4
|
Special versions of minCoeff and maxCoeff :
Typical use cases of all() and any():
Read-write access to a column or a row of a matrix (or array):
Read-write access to sub-vectors:
Default versions | Optimized versions when the size is known at compile time | |
---|---|---|
the first n coeffs | ||
the last n coeffs | ||
the n coeffs in the range [ pos : pos + n - 1] | ||
Read-write access to sub-matrices: | ||
(more) | (more) | the rows x cols sub-matrix starting from position ( i ,j ) |
the rows x cols sub-matrix taken in one of the four corners | ||
specialized versions of block() when the block fit two corners |
Vectors, rows, and/or columns of a matrix can be reversed (see DenseBase::reverse(), DenseBase::reverseInPlace(), VectorwiseOp::reverse()).
Vectors, matrices, rows, and/or columns can be replicated in any direction (see DenseBase::replicate(), VectorwiseOp::replicate())
(matrix world *)
Operation | Code |
---|---|
view a vector as a diagonal matrix | const DiagonalWrapper< const Derived > asDiagonal() const Definition: DiagonalMatrix.h:384 |
Declare a diagonal matrix | DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
diag1.diagonal() = vector;
|
Access the diagonal and super/sub diagonals of a matrix as a vector (read/write) | |
Optimized products and inverse | mat3 = scalar * diag1 * mat1;
|
TriangularView gives a view on a triangular part of a dense matrix and allows to perform optimized operations on it. The opposite triangular part is never referenced and can be used to store other information.
template
keyword if it is used on an object of a type that depends on a template parameter; see The template and typename keywords in C++ for details.Operation | Code |
---|---|
Reference to a triangular with optional unit or null diagonal (read/write): | m.triangularView<Xxx>()
TriangularViewReturnType< Mode >::Type triangularView() Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower |
Writing to a specific triangular part: (only the referenced triangular part is evaluated) | |
Conversion to a dense matrix setting the opposite triangular part to zero: | |
Products: | |
Solving linear equations: \( M_2 := L_1^{-1} M_2 \) \( M_3 := {L_1^*}^{-1} M_3 \) \( M_4 := M_4 U_1^{-1} \) | L1.triangularView<Eigen::UnitLower>().solveInPlace(M2)
L1.triangularView<Eigen::Lower>().adjoint().solveInPlace(M3)
U1.triangularView<Eigen::Upper>().solveInPlace<OnTheRight>(M4)
M1<< 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12;Map< MatrixXf > M2(M1.data(), 6, 2) RowMajorMatrixXf M3(M1) cout<< "Row major input:"<< endl<< M3<< "\n";Map< RowMajorMatrixXf, 0, Stride< Dynamic, 3 > > M4(M3.data(), M3.rows(),(M3.cols()+2)/3, Stride< Dynamic, 3 >(M3.outerStride(), 3)) |
Just as for triangular matrix, you can reference any triangular part of a square matrix to see it as a selfadjoint matrix and perform special and optimized operations. Again the opposite triangular part is never referenced and can be used to store other information.
template
keyword if it is used on an object of a type that depends on a template parameter; see The template and typename keywords in C++ for details.Operation | Code |
---|---|
Conversion to a dense matrix: | SelfAdjointViewReturnType< UpLo >::Type selfadjointView() |
Product with another general matrix or vector: | |
Rank 1 and rank K update: \( upper(M_1) \mathrel{{+}{=}} s_1 M_2 M_2^* \) \( lower(M_1) \mathbin{{-}{=}} M_2^* M_2 \) | |
Rank 2 update: ( \( M \mathrel{{+}{=}} s u v^* + s v u^* \)) | |
Solving linear equations: ( \( M_2 := M_1^{-1} M_2 \)) | // via a standard Cholesky factorization
// via a Cholesky factorization with pivoting
|