namespace Eigen { /** \page QuickRefPage Quick reference guide \b Table \b of \b contents - \ref QuickRef_Headers - \ref QuickRef_Types - \ref QuickRef_Map - \ref QuickRef_ArithmeticOperators - \ref QuickRef_Coeffwise - \ref QuickRef_Reductions - \ref QuickRef_Blocks - \ref QuickRef_DiagTriSymm \n
top \section QuickRef_Headers Modules and Header files
ModuleHeader fileContents
Core\code#include \endcodeMatrix and Array classes, basic linear algebra (including triangular and selfadjoint products), array manipulation
Geometry\code#include \endcodeTransformation, Translation, Scaling, 2D and 3D rotations (Quaternion, AngleAxis)
LU\code#include \endcodeInverse, determinant, LU decompositions (FullPivLU, PartialPivLU) with solver
Cholesky\code#include \endcodeLLT and LDLT Cholesky factorization with solver
SVD\code#include \endcodeSVD decomposition with solver (HouseholderSVD, JacobiSVD)
QR\code#include \endcodeQR decomposition with solver (HouseholderQR, ColPivHouseholerQR, FullPivHouseholderQR)
Eigenvalues\code#include \endcodeEigenvalue, eigenvector decompositions for selfadjoint and non selfadjoint real or complex matrices.
Sparse\code#include \endcodeSparse matrix storage and related basic linear algebra.
\code#include \endcodeIncludes Core, Geometry, LU, Cholesky, SVD, QR, and Eigenvalues
\code#include \endcodeIncludes Dense and Sparse
top \section QuickRef_Types Array, matrix and vector types \b 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: \code typedef Matrix MyMatrixType; typedef Array MyArrayType; \endcode \li \c Scalar is the scalar type of the coefficients (e.g., \c float, \c double, \c bool, \c int, etc.). \li \c RowsAtCompileTime and \c ColsAtCompileTime are the number of rows and columns of the matrix as known at compile-time or \c Dynamic. \li \c Options can be \c ColMajor or \c RowMajor, default is \c 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: \code Matrix // Dynamic number of columns (heap allocation) Matrix // Dynamic number of rows (heap allocation) Matrix // Fully dynamic, row major (heap allocation) Matrix // Fully fixed (static allocation) \endcode In most cases, you can simply use one of the convenience typedefs for \ref matrixtypedefs "matrices" and \ref arraytypedefs "arrays". Some examples:
\code Matrix <=> MatrixXf Matrix <=> VectorXd Matrix <=> RowVectorXi Matrix <=> Matrix3f Matrix <=> Vector4f \endcode\code Array <=> ArrayXXf Array <=> ArrayXd Array <=> RowArrayXi Array <=> Array33f Array <=> Array4f \endcode
Conversion between the matrix and array worlds: \code Array44f a1, a1; Matrix4f m1, m2; m1 = a1 * a2; // OK, coeffwise product a1 = m1 * m2; // OK, matrix product a2 = a1 + m1.array(); m2 = a1.matrix() + m1; ArrayWrapper m1a(m1); // m1a is an alias for m1.array(), they share the same coefficients MatrixWrapper a1m(a1); \endcode 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: \li \matrixworld linear algebra matrix and vector only \li \arrayworld array objects only \subsection QuickRef_Basics Basic matrix manipulation no-op if the new sizes match,\n otherwise data are lost \n \n resizing with data preservation
1D objects2D objectsNotes
Constructors \code Vector4d v4; Vector2f v1(x, y); Array3i v2(x, y, z); Vector4d v3(x, y, z, w); VectorXf v5; ArrayXf v6(size); \endcode\code Matrix4f m1; MatrixXf m5; MatrixXf m6(nb_rows, nb_columns); \endcode The coeffs are left uninitialized \n \n \n \n Empty object \n The coeffs are left uninitialized
Comma initializer \code Vector3f v1; v1 << x, y, z; ArrayXf v2(4); v2 << 1, 2, 3, 4; \endcode\code Matrix3f m1; m1 << 1, 2, 3, 4, 5, 6, 7, 8, 9; \endcode
Comma initializer (bis)
\include Tutorial_commainit_02.cpp output: \verbinclude Tutorial_commainit_02.out
Runtime info \code vector.size(); vector.innerStride(); vector.data(); \endcode\code matrix.rows(); matrix.cols(); matrix.innerSize(); matrix.outerSize(); matrix.innerStride(); matrix.outerStride(); matrix.data(); \endcode\n Inner/Outer* are storage order dependent
Compile-time info \code ObjectType::Scalar ObjectType::RowsAtCompileTime ObjectType::RealScalar ObjectType::ColsAtCompileTime ObjectType::Index ObjectType::SizeAtCompileTime \endcode
Resizing \code vector.resize(size); vector.resizeLike(other_vector); vector.conservativeResize(size); \endcode\code 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); \endcode
Coeff access with \n range checking \code vector(i) vector.x() vector[i] vector.y() vector.z() vector.w() \endcode\code matrix(i,j) \endcodeRange checking is disabled if \n NDEBUG or EIGEN_NO_DEBUG is defined
Coeff access without \n range checking \code vector.coeff(i) vector.coeffRef(i) \endcode\code matrix.coeff(i,j) matrix.coeffRef(i,j) \endcode
Assignment/copy \code object = expression; object_of_float = expression_of_double.cast(); \endcodethe destination is automatically resized (if possible)
\subsection QuickRef_PredefMat Predefined Matrices
Fixed-size matrix or vector Dynamic-size matrix Dynamic-size vector
\code typedef {Matrix3f|Array33f} FixedXD; FixedXD x; x = FixedXD::Zero(); x = FixedXD::Ones(); x = FixedXD::Constant(value); x = FixedXD::Random(); x.setZero(); x.setOnes(); x.setConstant(value); x.setRandom(); \endcode \code typedef {MatrixXf|ArrayXXf} Dynamic2D; Dynamic2D x; x = Dynamic2D::Zero(rows, cols); x = Dynamic2D::Ones(rows, cols); x = Dynamic2D::Constant(rows, cols, value); x = Dynamic2D::Random(rows, cols); x.setZero(rows, cols); x.setOnes(rows, cols); x.setConstant(rows, cols, value); x.setRandom(rows, cols); \endcode \code typedef {VectorXf|ArrayXf} Dynamic1D; Dynamic1D x; x = Dynamic1D::Zero(size); x = Dynamic1D::Ones(size); x = Dynamic1D::Constant(size, value); x = Dynamic1D::Random(size); x.setZero(size); x.setOnes(size); x.setConstant(size, value); x.setRandom(size); \endcode
Identity and \link MatrixBase::Unit basis vectors \endlink \matrixworld
\code x = FixedXD::Identity(); x.setIdentity(); Vector3f::UnitX() // 1 0 0 Vector3f::UnitY() // 0 1 0 Vector3f::UnitZ() // 0 0 1 \endcode \code x = Dynamic2D::Identity(rows, cols); x.setIdentity(rows, cols); N/A \endcode \code N/A VectorXf::Unit(size,i) VectorXf::Unit(4,1) == Vector4f(0,1,0,0) == Vector4f::UnitY() \endcode
\subsection QuickRef_Map Map
Contiguous memory \code float data[] = {1,2,3,4}; Map v1(data); // uses v1 as a Vector3f object Map v2(data,3); // uses v2 as a ArrayXf object Map m1(data); // uses m1 as a Array22f object Map m2(data,2,2); // uses m2 as a MatrixXf object \endcode
Typical usage of strides \code float data[] = {1,2,3,4,5,6,7,8,9}; Map > v1(data,3); // == [1,3,5] Map > v2(data,3,InnerStride<>(3)); // == [1,4,7] Map > m1(data,2,3,OuterStride<>(3)); // == |1,4,7| Map > m2(data,2,3); // |2,5,8| \endcode
top \section QuickRef_ArithmeticOperators Arithmetic Operators
add/subtract\code mat3 = mat1 + mat2; mat3 += mat1; mat3 = mat1 - mat2; mat3 -= mat1;\endcode
scalar product\code mat3 = mat1 * s1; mat3 = s1 * mat1; mat3 *= s1; mat3 = mat1 / s1; mat3 /= s1;\endcode
matrix/vector product \matrixworld\code col2 = mat1 * col1; row2 = row1 * mat1; row1 *= mat1; mat3 = mat1 * mat2; mat3 *= mat1; \endcode
transpose et adjoint \matrixworld\code mat1 = mat2.transpose(); mat1.transposeInPlace(); mat1 = mat2.adjoint(); mat1.adjointInPlace(); \endcode
\link MatrixBase::dot() dot \endlink \& inner products \matrixworld\code scalar = col1.adjoint() * col2; scalar = (col1.adjoint() * col2).value(); scalar = vec1.dot(vec2);\endcode
outer product \matrixworld\code mat = col1 * col2.transpose();\endcode
\link MatrixBase::norm() norm \endlink and \link MatrixBase::normalized() normalization \endlink \matrixworld\code scalar = vec1.norm(); scalar = vec1.squaredNorm() vec2 = vec1.normalized(); vec1.normalize(); // inplace \endcode
\link MatrixBase::cross() cross product \endlink \matrixworld\code #include vec3 = vec1.cross(vec2);\endcode
top \section QuickRef_Coeffwise Coefficient-wise \& Array operators Coefficient-wise operators for matrices and vectors:
Matrix API \matrixworldVia Array conversions
\code mat1.cwiseMin(mat2) mat1.cwiseMax(mat2) mat1.cwiseAbs2() mat1.cwiseAbs() mat1.cwiseSqrt() mat1.cwiseProduct(mat2) mat1.cwiseQuotient(mat2)\endcode \code mat1.array().min(mat2.array()) mat1.array().max(mat2.array()) mat1.array().abs2() mat1.array().abs() mat1.array().sqrt() mat1.array() * mat2.array() mat1.array() / mat2.array() \endcode
Array operators:\arrayworld
Arithmetic operators\code array1 * array2 array1 / array2 array1 *= array2 array1 /= array2 array1 + scalar array1 - scalar array1 += scalar array1 -= scalar \endcode
Comparisons\code array1 < array2 array1 > array2 array1 < scalar array1 > scalar array1 <= array2 array1 >= array2 array1 <= scalar array1 >= scalar array1 == array2 array1 != array2 array1 == scalar array1 != scalar \endcode
Special functions \n and STL variants\code array1.min(array2) std::min(array1,array2) array1.max(array2) std::max(array1,array2) array1.abs2() array1.abs() std::abs(array1) array1.sqrt() std::sqrt(array1) array1.log() std::log(array1) array1.exp() std::exp(array1) array1.pow(exponent) std::pow(array1,exponent) array1.square() array1.cube() array1.inverse() array1.sin() std::sin(array1) array1.cos() std::cos(array1) array1.tan() std::tan(array1) \endcode
top \section QuickRef_Reductions Reductions Eigen provides several reduction methods such as: \link DenseBase::minCoeff() minCoeff() \endlink, \link DenseBase::maxCoeff() maxCoeff() \endlink, \link DenseBase::sum() sum() \endlink, \link MatrixBase::trace() trace() \endlink \matrixworld, \link MatrixBase::norm() norm() \endlink \matrixworld, \link MatrixBase::squaredNorm() squaredNorm() \endlink \matrixworld, \link DenseBase::all() all() \endlink \redstar,and \link DenseBase::any() any() \endlink \redstar. All reduction operations can be done matrix-wise, \link DenseBase::colwise() column-wise \endlink \redstar or \link DenseBase::rowwise() row-wise \endlink \redstar. Usage example:
\code 5 3 1 mat = 2 7 8 9 4 6 \endcode \code mat.minCoeff(); \endcode\code 1 \endcode
\code mat.colwise().minCoeff(); \endcode\code 2 3 1 \endcode
\code mat.rowwise().minCoeff(); \endcode\code 1 2 4 \endcode
Also note that maxCoeff and minCoeff can takes optional arguments returning the coordinates of the respective min/max coeff: \link DenseBase::maxCoeff(int*,int*) const maxCoeff(int* i, int* j) \endlink, \link DenseBase::minCoeff(int*,int*) const minCoeff(int* i, int* j) \endlink. \b Side \b note: The all() and any() functions are especially useful in combination with coeff-wise comparison operators. top\section QuickRef_Blocks Matrix blocks Read-write access to a \link DenseBase::col(int) column \endlink or a \link DenseBase::row(int) row \endlink of a matrix (or array): \code mat1.row(i) = mat2.col(j); mat1.col(j1).swap(mat1.col(j2)); \endcode Read-write access to sub-vectors:
Default versions Optimized versions when the size \n is known at compile time
\code vec1.head(n)\endcode\code vec1.head()\endcodethe first \c n coeffs
\code vec1.tail(n)\endcode\code vec1.tail()\endcodethe last \c n coeffs
\code vec1.segment(pos,n)\endcode\code vec1.segment(pos)\endcode the \c size coeffs in \n the range [\c pos : \c pos + \c n [
Read-write access to sub-matrices:
\code mat1.block(i,j,rows,cols)\endcode \link DenseBase::block(int,int,int,int) (more) \endlink \code mat1.block(i,j)\endcode \link DenseBase::block(int,int) (more) \endlink the \c rows x \c cols sub-matrix \n starting from position (\c i,\c j)
\code mat1.topLeftCorner(rows,cols) mat1.topRightCorner(rows,cols) mat1.bottomLeftCorner(rows,cols) mat1.bottomRightCorner(rows,cols)\endcode \code mat1.topLeftCorner() mat1.topRightCorner() mat1.bottomLeftCorner() mat1.bottomRightCorner()\endcode the \c rows x \c cols sub-matrix \n taken in one of the four corners
\code mat1.topRows(rows) mat1.bottomRows(rows) mat1.leftCols(cols) mat1.rightCols(cols)\endcode \code mat1.topRows() mat1.bottomRows() mat1.leftCols() mat1.rightCols()\endcode specialized versions of block() when the block fit two corners
top\section QuickRef_DiagTriSymm Diagonal, Triangular, and Self-adjoint matrices (matrix world \matrixworld) \subsection QuickRef_Diagonal Diagonal matrices
\link MatrixBase::asDiagonal() make a diagonal matrix \endlink \n from a vector \code mat1 = vec1.asDiagonal();\endcode
Declare a diagonal matrix\code DiagonalMatrix diag1(size); diag1.diagonal() = vector;\endcode
Access \link MatrixBase::diagonal() the diagonal and super/sub diagonals of a matrix \endlink as a vector (read/write) \code vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal \endcode
Optimized products and inverse \code mat3 = scalar * diag1 * mat1; mat3 += scalar * mat1 * vec1.asDiagonal(); mat3 = vec1.asDiagonal().inverse() * mat1 mat3 = mat1 * diag1.inverse() \endcode
\subsection QuickRef_TriangularView Triangular views TriangularView allows to get views on a triangular part of a dense matrix and perform optimized operations on it. The opposite triangular is never referenced and can be used to store other information.
Reference a read/write triangular part of a given \n matrix (or expression) m with optional unit diagonal: \code m.triangularView() \endcode \n \c Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
Writing to a specific triangular part:\n (only the referenced triangular part is evaluated) \code m1.triangularView() = m2 + m3 \endcode
Conversion to a dense matrix setting the opposite triangular part to zero: \code m2 = m1.triangularView()\endcode
Products: \code m3 += s1 * m1.adjoint().triangularView() * m2 m3 -= s1 * m2.conjugate() * m1.adjoint().triangularView() \endcode
Solving linear equations:\n(\f$ m_2 := m_1^{-1} m_2 \f$) \code m1.triangularView().solveInPlace(m2) m1.adjoint().triangularView().solveInPlace(m2)\endcode
\subsection QuickRef_SelfadjointMatrix Symmetric/selfadjoint views Just as for triangular matrix, you can reference any triangular part of a square matrix to see it a selfadjoint matrix and perform special and optimized operations. Again the opposite triangular is never referenced and can be used to store other information.
Conversion to a dense matrix: \code m2 = m.selfadjointView();\endcode
Product with another general matrix or vector: \code m3 = s1 * m1.conjugate().selfadjointView() * m3; m3 -= s1 * m3.adjoint() * m1.selfadjointView();\endcode
Rank 1 and rank K update: \code // fast version of m1 += s1 * m2 * m2.adjoint(): m1.selfadjointView().rankUpdate(m2,s1); // fast version of m1 -= m2.adjoint() * m2: m1.selfadjointView().rankUpdate(m2.adjoint(),-1); \endcode
Rank 2 update: (\f$ m += s u v^* + s v u^* \f$) \code m.selfadjointView().rankUpdate(u,v,s); \endcode
Solving linear equations:\n(\f$ m_2 := m_1^{-1} m_2 \f$) \code // via a standard Cholesky factorization m1.selfadjointView().llt().solveInPlace(m2); // via a Cholesky factorization with pivoting m1.selfadjointView().ldlt().solveInPlace(m2); \endcode
*/ /*
\link MatrixBase::asDiagonal() make a diagonal matrix \endlink \n from a vector \code mat1 = vec1.asDiagonal();\endcode
Declare a diagonal matrix\code DiagonalMatrix diag1(size); diag1.diagonal() = vector;\endcode
Access \link MatrixBase::diagonal() the diagonal and super/sub diagonals of a matrix \endlink as a vector (read/write) \code vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal \endcode
View on a triangular part of a matrix (read/write) \code mat2 = mat1.triangularView(); // Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower mat1.triangularView() = mat2 + mat3; // only the upper part is evaluated and referenced \endcode
View a triangular part as a symmetric/self-adjoint matrix (read/write) \code mat2 = mat1.selfadjointView(); // Xxx = Upper or Lower mat1.selfadjointView() = mat2 + mat2.adjoint(); // evaluated and write to the upper triangular part only \endcode
Optimized products: \code mat3 += scalar * vec1.asDiagonal() * mat1 mat3 += scalar * mat1 * vec1.asDiagonal() mat3.noalias() += scalar * mat1.triangularView() * mat2 mat3.noalias() += scalar * mat2 * mat1.triangularView() mat3.noalias() += scalar * mat1.selfadjointView() * mat2 mat3.noalias() += scalar * mat2 * mat1.selfadjointView() mat1.selfadjointView().rankUpdate(mat2); mat1.selfadjointView().rankUpdate(mat2.adjoint(), scalar); \endcode Inverse products: (all are optimized) \code mat3 = vec1.asDiagonal().inverse() * mat1 mat3 = mat1 * diag1.inverse() mat1.triangularView().solveInPlace(mat2) mat1.triangularView().solveInPlace(mat2) mat2 = mat1.selfadjointView().llt().solve(mat2) \endcode */ }