namespace Eigen { /** \eigenManualPage QuickRefPage Quick reference guide \eigenAutoToc
top \section QuickRef_Headers Modules and Header files 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 \c %Dense and \c Eigen header files are provided to conveniently gain access to several modules at once.
ModuleHeader fileContents
\link Core_Module Core \endlink\code#include \endcodeMatrix and Array classes, basic linear algebra (including triangular and selfadjoint products), array manipulation
\link Geometry_Module Geometry \endlink\code#include \endcodeTransform, Translation, Scaling, Rotation2D and 3D rotations (Quaternion, AngleAxis)
\link LU_Module LU \endlink\code#include \endcodeInverse, determinant, LU decompositions with solver (FullPivLU, PartialPivLU)
\link Cholesky_Module Cholesky \endlink\code#include \endcodeLLT and LDLT Cholesky factorization with solver
\link Householder_Module Householder \endlink\code#include \endcodeHouseholder transformations; this module is used by several linear algebra modules
\link SVD_Module SVD \endlink\code#include \endcodeSVD decompositions with least-squares solver (JacobiSVD, BDCSVD)
\link QR_Module QR \endlink\code#include \endcodeQR decomposition with solver (HouseholderQR, ColPivHouseholderQR, FullPivHouseholderQR)
\link Eigenvalues_Module Eigenvalues \endlink\code#include \endcodeEigenvalue, eigenvector decompositions (EigenSolver, SelfAdjointEigenSolver, ComplexEigenSolver)
\link Sparse_Module Sparse \endlink\code#include \endcode%Sparse matrix storage and related basic linear algebra (SparseMatrix, SparseVector) \n (see \ref SparseQuickRefPage for details on sparse modules)
\code#include \endcodeIncludes Core, Geometry, LU, Cholesky, SVD, QR, and Eigenvalues header files
\code#include \endcodeIncludes %Dense and %Sparse header files (the whole Eigen library)
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 (usually allocated on stack) \endcode In most cases, you can simply use one of the convenience typedefs for \ref matrixtypedefs "matrices" and \ref arraytypedefs "arrays". Some examples:
MatricesArrays
\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, a2; Matrix4f m1, m2; m1 = a1 * a2; // coeffwise product, implicit conversion from array to matrix. a1 = m1 * m2; // matrix product, implicit conversion from matrix to array. a2 = a1 + m1.array(); // mixing array and matrix is forbidden m2 = a1.matrix() + m1; // and explicit conversion is required. 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
1D objects2D objectsNotes
Constructors \code Vector4d v4; Vector2f v1(x, y); Array3i v2(x, y, z); Vector4d v3(x, y, z, w); VectorXf v5; // empty object ArrayXf v6(size); \endcode\code Matrix4f m1; MatrixXf m5; // empty object MatrixXf m6(nb_rows, nb_columns); \endcode By default, the coefficients \n 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(); \endcodeInner/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); \endcodeno-op if the new sizes match,
otherwise data are lost

resizing with data preservation
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 = FixedXD::LinSpaced(size, low, high); x.setZero(); x.setOnes(); x.setConstant(value); x.setRandom(); x.setLinSpaced(size, low, high); \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); N/A x.setZero(rows, cols); x.setOnes(rows, cols); x.setConstant(rows, cols, value); x.setRandom(rows, cols); N/A \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 = Dynamic1D::LinSpaced(size, low, high); x.setZero(size); x.setOnes(size); x.setConstant(size, value); x.setRandom(size); x.setLinSpaced(size, low, high); \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 Vector4f::Unit(i) x.setUnit(i); \endcode \code x = Dynamic2D::Identity(rows, cols); x.setIdentity(rows, cols); N/A \endcode \code N/A VectorXf::Unit(size,i) x.setUnit(size,i); VectorXf::Unit(4,1) == Vector4f(0,1,0,0) == Vector4f::UnitY() \endcode
Note that it is allowed to call any of the \c set* functions to a dynamic-sized vector or matrix without passing new sizes. For instance: \code MatrixXi M(3,3); M.setIdentity(); \endcode \subsection QuickRef_Map Mapping external arrays
Contiguous \n 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 \n 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 > m2(data,2,3); // both lines |1,4,7| Map > m1(data,2,3,OuterStride<>(3)); // are equal to: |2,5,8| \endcode
top \section QuickRef_ArithmeticOperators Arithmetic Operators
add \n subtract\code mat3 = mat1 + mat2; mat3 += mat1; mat3 = mat1 - mat2; mat3 -= mat1;\endcode
scalar product\code mat3 = mat1 * s1; mat3 *= s1; mat3 = s1 * mat1; mat3 = mat1 / s1; mat3 /= s1;\endcode
matrix/vector \n products \matrixworld\code col2 = mat1 * col1; row2 = row1 * mat1; row1 *= mat1; mat3 = mat1 * mat2; mat3 *= mat1; \endcode
transposition \n adjoint \matrixworld\code mat1 = mat2.transpose(); mat1.transposeInPlace(); mat1 = mat2.adjoint(); mat1.adjointInPlace(); \endcode
\link MatrixBase::dot dot \endlink product \n inner product \matrixworld\code scalar = vec1.dot(vec2); scalar = col1.adjoint() * col2; scalar = (col1.adjoint() * col2).value();\endcode
outer product \matrixworld\code mat = col1 * col2.transpose();\endcode
\link MatrixBase::norm() norm \endlink \n \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 In addition to the aforementioned operators, Eigen supports numerous coefficient-wise operator and functions. Most of them unambiguously makes sense in array-world\arrayworld. The following operators are readily available for arrays, or available through .array() for vectors and matrices:
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 array1.min(array2) array1.max(array2) array1.min(scalar) array1.max(scalar) \endcode
Trigo, power, and \n misc functions \n and the STL-like variants\code 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) \endcode
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 APISTL-like APIs\arrayworld Comments
\code mat1.real() mat1.imag() mat1.conjugate() \endcode \code real(array1) imag(array1) conj(array1) \endcode \code // read-write, no-op for real expressions // read-only for real, read-write for complexes // no-op for real expressions \endcode
Some coefficient-wise operators are readily available for for matrices and vectors through the following cwise* methods:
Matrix API \matrixworldVia Array conversions
\code mat1.cwiseMin(mat2) mat1.cwiseMin(scalar) mat1.cwiseMax(mat2) mat1.cwiseMax(scalar) mat1.cwiseAbs2() mat1.cwiseAbs() mat1.cwiseSqrt() mat1.cwiseInverse() mat1.cwiseProduct(mat2) mat1.cwiseQuotient(mat2) mat1.cwiseEqual(mat2) mat1.cwiseEqual(scalar) mat1.cwiseNotEqual(mat2) \endcode \code mat1.array().min(mat2.array()) mat1.array().min(scalar) mat1.array().max(mat2.array()) mat1.array().max(scalar) mat1.array().abs2() mat1.array().abs() mat1.array().sqrt() mat1.array().inverse() mat1.array() * mat2.array() mat1.array() / mat2.array() mat1.array() == mat2.array() mat1.array() == scalar mat1.array() != mat2.array() \endcode
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 \c foo using DenseBase::unaryExpr together with std::ptr_fun (c++03), std::ref (c++11), or lambdas (c++11): \code mat1.unaryExpr(std::ptr_fun(foo)); mat1.unaryExpr(std::ref(foo)); mat1.unaryExpr([](double x) { return foo(x); }); \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 DenseBase::prod() prod() \endlink, \link MatrixBase::trace() trace() \endlink \matrixworld, \link MatrixBase::norm() norm() \endlink \matrixworld, \link MatrixBase::squaredNorm() squaredNorm() \endlink \matrixworld, \link DenseBase::all() all() \endlink, and \link DenseBase::any() any() \endlink. All reduction operations can be done matrix-wise, \link DenseBase::colwise() column-wise \endlink or \link DenseBase::rowwise() row-wise \endlink. 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
Special versions of \link DenseBase::minCoeff(IndexType*,IndexType*) const minCoeff \endlink and \link DenseBase::maxCoeff(IndexType*,IndexType*) const maxCoeff \endlink: \code int i, j; s = vector.minCoeff(&i); // s == vector[i] s = matrix.maxCoeff(&i, &j); // s == matrix(i,j) \endcode Typical use cases of all() and any(): \code if((array1 > 0).all()) ... // if all coefficients of array1 are greater than 0 ... if((array1 < array2).any()) ... // if there exist a pair i,j such that array1(i,j) < array2(i,j) ... \endcode top\section QuickRef_Blocks Sub-matrices Read-write access to a \link DenseBase::col(Index) column \endlink or a \link DenseBase::row(Index) 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 n coeffs in the \n range [\c pos : \c pos + \c n - 1]
Read-write access to sub-matrices:
\code mat1.block(i,j,rows,cols)\endcode \link DenseBase::block(Index,Index,Index,Index) (more) \endlink \code mat1.block(i,j)\endcode \link DenseBase::block(Index,Index) (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() \n when the block fit two corners
top\section QuickRef_Misc Miscellaneous operations \subsection QuickRef_Reverse Reverse Vectors, rows, and/or columns of a matrix can be reversed (see DenseBase::reverse(), DenseBase::reverseInPlace(), VectorwiseOp::reverse()). \code vec.reverse() mat.colwise().reverse() mat.rowwise().reverse() vec.reverseInPlace() \endcode \subsection QuickRef_Replicate Replicate Vectors, matrices, rows, and/or columns can be replicated in any direction (see DenseBase::replicate(), VectorwiseOp::replicate()) \code vec.replicate(times) vec.replicate mat.replicate(vertical_times, horizontal_times) mat.replicate() mat.colwise().replicate(vertical_times, horizontal_times) mat.colwise().replicate() mat.rowwise().replicate(vertical_times, horizontal_times) mat.rowwise().replicate() \endcode top\section QuickRef_DiagTriSymm Diagonal, Triangular, and Self-adjoint matrices (matrix world \matrixworld) \subsection QuickRef_Diagonal Diagonal matrices
OperationCode
view a vector \link MatrixBase::asDiagonal() as a diagonal matrix \endlink \n \code mat1 = vec1.asDiagonal();\endcode
Declare a diagonal matrix\code DiagonalMatrix diag1(size); diag1.diagonal() = vector;\endcode
Access the \link MatrixBase::diagonal() diagonal \endlink and \link MatrixBase::diagonal(Index) super/sub diagonals \endlink of a matrix 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 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. \note The .triangularView() template member function requires the \c template keyword if it is used on an object of a type that depends on a template parameter; see \ref TopicTemplateKeyword for details.
OperationCode
Reference to a triangular with optional \n unit or null diagonal (read/write): \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 := L_1^{-1} M_2 \f$ \n \f$ M_3 := {L_1^*}^{-1} M_3 \f$ \n \f$ M_4 := M_4 U_1^{-1} \f$ \n \code L1.triangularView().solveInPlace(M2) L1.triangularView().adjoint().solveInPlace(M3) U1.triangularView().solveInPlace(M4)\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 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. \note The .selfadjointView() template member function requires the \c template keyword if it is used on an object of a type that depends on a template parameter; see \ref TopicTemplateKeyword for details.
OperationCode
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: \n \f$ upper(M_1) \mathrel{{+}{=}} s_1 M_2 M_2^* \f$ \n \f$ lower(M_1) \mathbin{{-}{=}} M_2^* M_2 \f$ \n \code M1.selfadjointView().rankUpdate(M2,s1); M1.selfadjointView().rankUpdate(M2.adjoint(),-1); \endcode
Rank 2 update: (\f$ M \mathrel{{+}{=}} 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 m2 = m1.selfadjointView().llt().solve(m2); // via a Cholesky factorization with pivoting m2 = m1.selfadjointView().ldlt().solve(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 */ }